xref: /petsc/src/snes/impls/vi/vi.c (revision 69c03620c67e0552d305d2275534bbb82c27e56b)
12b4ed282SShri Abhyankar #define PETSCSNES_DLL
22b4ed282SShri Abhyankar 
32b4ed282SShri Abhyankar #include "../src/snes/impls/vi/viimpl.h"
42b4ed282SShri Abhyankar 
52b4ed282SShri Abhyankar /*
62b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
72b4ed282SShri Abhyankar     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
82b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
92b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
102b4ed282SShri Abhyankar */
112b4ed282SShri Abhyankar #undef __FUNCT__
122b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
13ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
142b4ed282SShri Abhyankar {
152b4ed282SShri Abhyankar   PetscReal      a1;
162b4ed282SShri Abhyankar   PetscErrorCode ierr;
17ace3abfcSBarry Smith   PetscBool     hastranspose;
182b4ed282SShri Abhyankar 
192b4ed282SShri Abhyankar   PetscFunctionBegin;
202b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
212b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
222b4ed282SShri Abhyankar   if (hastranspose) {
232b4ed282SShri Abhyankar     /* Compute || J^T F|| */
242b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
252b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
262b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
272b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
282b4ed282SShri Abhyankar   } else {
292b4ed282SShri Abhyankar     Vec         work;
302b4ed282SShri Abhyankar     PetscScalar result;
312b4ed282SShri Abhyankar     PetscReal   wnorm;
322b4ed282SShri Abhyankar 
332b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
342b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
352b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
362b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
372b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
382b4ed282SShri Abhyankar     ierr = VecDestroy(work);CHKERRQ(ierr);
392b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
402b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
412b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
422b4ed282SShri Abhyankar   }
432b4ed282SShri Abhyankar   PetscFunctionReturn(0);
442b4ed282SShri Abhyankar }
452b4ed282SShri Abhyankar 
462b4ed282SShri Abhyankar /*
472b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
482b4ed282SShri Abhyankar */
492b4ed282SShri Abhyankar #undef __FUNCT__
502b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
512b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
522b4ed282SShri Abhyankar {
532b4ed282SShri Abhyankar   PetscReal      a1,a2;
542b4ed282SShri Abhyankar   PetscErrorCode ierr;
55ace3abfcSBarry Smith   PetscBool     hastranspose;
562b4ed282SShri Abhyankar 
572b4ed282SShri Abhyankar   PetscFunctionBegin;
582b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
592b4ed282SShri Abhyankar   if (hastranspose) {
602b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
612b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
622b4ed282SShri Abhyankar 
632b4ed282SShri Abhyankar     /* Compute || J^T W|| */
642b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
652b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
662b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
672b4ed282SShri Abhyankar     if (a1 != 0.0) {
682b4ed282SShri Abhyankar       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
692b4ed282SShri Abhyankar     }
702b4ed282SShri Abhyankar   }
712b4ed282SShri Abhyankar   PetscFunctionReturn(0);
722b4ed282SShri Abhyankar }
732b4ed282SShri Abhyankar 
742b4ed282SShri Abhyankar /*
752b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
762b4ed282SShri Abhyankar 
772b4ed282SShri Abhyankar   Notes:
782b4ed282SShri Abhyankar   The convergence criterion currently implemented is
792b4ed282SShri Abhyankar   merit < abstol
802b4ed282SShri Abhyankar   merit < rtol*merit_initial
812b4ed282SShri Abhyankar */
822b4ed282SShri Abhyankar #undef __FUNCT__
832b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
842b4ed282SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy)
852b4ed282SShri Abhyankar {
862b4ed282SShri Abhyankar   PetscErrorCode ierr;
872b4ed282SShri Abhyankar 
882b4ed282SShri Abhyankar   PetscFunctionBegin;
892b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
902b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
912b4ed282SShri Abhyankar 
922b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
932b4ed282SShri Abhyankar 
942b4ed282SShri Abhyankar   if (!it) {
952b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
962b4ed282SShri Abhyankar     snes->ttol = merit*snes->rtol;
972b4ed282SShri Abhyankar   }
982b4ed282SShri Abhyankar   if (merit != merit) {
992b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
1002b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1012b4ed282SShri Abhyankar   } else if (merit < snes->abstol) {
1022b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr);
1032b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
1042b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
1052b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
1062b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
1072b4ed282SShri Abhyankar   }
1082b4ed282SShri Abhyankar 
1092b4ed282SShri Abhyankar   if (it && !*reason) {
1102b4ed282SShri Abhyankar     if (merit < snes->ttol) {
1112b4ed282SShri Abhyankar       ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr);
1122b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
1132b4ed282SShri Abhyankar     }
1142b4ed282SShri Abhyankar   }
1152b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1162b4ed282SShri Abhyankar }
1172b4ed282SShri Abhyankar 
1182b4ed282SShri Abhyankar /*
1192b4ed282SShri Abhyankar   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
1202b4ed282SShri Abhyankar 
1212b4ed282SShri Abhyankar   Input Parameter:
1222b4ed282SShri Abhyankar . phi - the semismooth function
1232b4ed282SShri Abhyankar 
1242b4ed282SShri Abhyankar   Output Parameter:
1252b4ed282SShri Abhyankar . merit - the merit function
1262b4ed282SShri Abhyankar . phinorm - ||phi||
1272b4ed282SShri Abhyankar 
1282b4ed282SShri Abhyankar   Notes:
1292b4ed282SShri Abhyankar   The merit function for the mixed complementarity problem is defined as
1302b4ed282SShri Abhyankar      merit = 0.5*phi^T*phi
1312b4ed282SShri Abhyankar */
1322b4ed282SShri Abhyankar #undef __FUNCT__
1332b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction"
1342b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
1352b4ed282SShri Abhyankar {
1362b4ed282SShri Abhyankar   PetscErrorCode ierr;
1372b4ed282SShri Abhyankar 
1382b4ed282SShri Abhyankar   PetscFunctionBegin;
1392b4ed282SShri Abhyankar   ierr = VecNormBegin(phi,NORM_2,phinorm);
1402b4ed282SShri Abhyankar   ierr = VecNormEnd(phi,NORM_2,phinorm);
1412b4ed282SShri Abhyankar 
1422b4ed282SShri Abhyankar   *merit = 0.5*(*phinorm)*(*phinorm);
1432b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1442b4ed282SShri Abhyankar }
1452b4ed282SShri Abhyankar 
1462b4ed282SShri Abhyankar /*
1472b4ed282SShri Abhyankar   ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation.
1482b4ed282SShri Abhyankar 
1492b4ed282SShri Abhyankar   Notes:
1502b4ed282SShri Abhyankar   The Fischer-Burmeister function is defined as
1512b4ed282SShri Abhyankar        ff(a,b) := sqrt(a*a + b*b) - a - b
1522b4ed282SShri Abhyankar   and is used reformulate a complementarity equation  as a semismooth equation.
1532b4ed282SShri Abhyankar */
1542b4ed282SShri Abhyankar 
1552b4ed282SShri Abhyankar #undef __FUNCT__
1562b4ed282SShri Abhyankar #define __FUNCT__ "ComputeFischerFunction"
1572b4ed282SShri Abhyankar static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff)
1582b4ed282SShri Abhyankar {
1592b4ed282SShri Abhyankar   PetscFunctionBegin;
1602b4ed282SShri Abhyankar   *ff = sqrt(a*a + b*b) - a - b;
1612b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1622b4ed282SShri Abhyankar }
1632b4ed282SShri Abhyankar 
1642b4ed282SShri Abhyankar /*
1651f79c099SShri Abhyankar    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
1662b4ed282SShri Abhyankar 
1672b4ed282SShri Abhyankar    Input Parameters:
1682b4ed282SShri Abhyankar .  snes - the SNES context
1692b4ed282SShri Abhyankar .  x - current iterate
1702b4ed282SShri Abhyankar .  functx - user defined function context
1712b4ed282SShri Abhyankar 
1722b4ed282SShri Abhyankar    Output Parameters:
1732b4ed282SShri Abhyankar .  phi - Semismooth function
1742b4ed282SShri Abhyankar 
1752b4ed282SShri Abhyankar    Notes:
1762b4ed282SShri Abhyankar    The result of this function is done by cases:
1772b4ed282SShri Abhyankar +  l[i] == -infinity, u[i] == infinity  -- phi[i] = -f[i]
1782b4ed282SShri Abhyankar .  l[i] == -infinity, u[i] finite       -- phi[i] = ss(u[i]-x[i], -f[i])
1792b4ed282SShri Abhyankar .  l[i] finite,       u[i] == infinity  -- phi[i] = ss(x[i]-l[i],  f[i])
1802b4ed282SShri Abhyankar .  l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u]))
1812b4ed282SShri Abhyankar -  otherwise l[i] == u[i] -- phi[i] = l[i] - x[i]
1822b4ed282SShri Abhyankar    ss is the semismoothing function used to reformulate the nonlinear equation in complementarity
1832b4ed282SShri Abhyankar    form to semismooth form
1842b4ed282SShri Abhyankar 
1852b4ed282SShri Abhyankar */
1862b4ed282SShri Abhyankar #undef __FUNCT__
1871f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction"
1881f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
1892b4ed282SShri Abhyankar {
1902b4ed282SShri Abhyankar   PetscErrorCode  ierr;
1912b4ed282SShri Abhyankar   SNES_VI       *vi = (SNES_VI*)snes->data;
1922b4ed282SShri Abhyankar   Vec             Xl = vi->xl,Xu = vi->xu,F = snes->vec_func;
1932b4ed282SShri Abhyankar   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u,t;
1942b4ed282SShri Abhyankar   PetscInt        i,nlocal;
1952b4ed282SShri Abhyankar 
1962b4ed282SShri Abhyankar   PetscFunctionBegin;
1972b4ed282SShri Abhyankar   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
1982b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
1992b4ed282SShri Abhyankar 
2002b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
2012b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
2022b4ed282SShri Abhyankar   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
2032b4ed282SShri Abhyankar   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
2042b4ed282SShri Abhyankar   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
2052b4ed282SShri Abhyankar 
2062b4ed282SShri Abhyankar   for (i=0;i < nlocal;i++) {
2072b4ed282SShri Abhyankar     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {
2082b4ed282SShri Abhyankar       phi_arr[i] = -f_arr[i];
2092b4ed282SShri Abhyankar     }
2102b4ed282SShri Abhyankar     else if (l[i] <= PETSC_VI_NINF) {
2112b4ed282SShri Abhyankar       t = u[i] - x_arr[i];
2122b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr);
2132b4ed282SShri Abhyankar       phi_arr[i] = -phi_arr[i];
2142b4ed282SShri Abhyankar     }
2152b4ed282SShri Abhyankar     else if (u[i] >= PETSC_VI_INF) {
2162b4ed282SShri Abhyankar       t = x_arr[i] - l[i];
2172b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr);
2182b4ed282SShri Abhyankar     }
2192b4ed282SShri Abhyankar     else if (l[i] == u[i]) {
2202b4ed282SShri Abhyankar       phi_arr[i] = l[i] - x_arr[i];
2212b4ed282SShri Abhyankar     }
2222b4ed282SShri Abhyankar     else {
2232b4ed282SShri Abhyankar       t = u[i] - x_arr[i];
2242b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);
2252b4ed282SShri Abhyankar       t = x_arr[i] - l[i];
2262b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]);
2272b4ed282SShri Abhyankar     }
2282b4ed282SShri Abhyankar   }
2292b4ed282SShri Abhyankar 
2302b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
2312b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
2322b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
2332b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
2342b4ed282SShri Abhyankar   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
2352b4ed282SShri Abhyankar 
2362b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2372b4ed282SShri Abhyankar }
2382b4ed282SShri Abhyankar 
2392b4ed282SShri Abhyankar /*
2401f79c099SShri Abhyankar    SNESVIComputeJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems.
2412b4ed282SShri Abhyankar 
2422b4ed282SShri Abhyankar    Input Parameters:
2432b4ed282SShri Abhyankar .  snes     - the SNES context
2442b4ed282SShri Abhyankar .  X        - the current iterate
2452b4ed282SShri Abhyankar .  vec_func - nonlinear function evaluated at x
2462b4ed282SShri Abhyankar 
2472b4ed282SShri Abhyankar    Output Parameters:
2482b4ed282SShri Abhyankar .  jac      - semismooth jacobian
2492b4ed282SShri Abhyankar .  jac_pre  - optional preconditioning matrix
2502b4ed282SShri Abhyankar .  flag     - flag passed on by SNESComputeJacobian.
2512b4ed282SShri Abhyankar .  jacctx   - user provided jacobian context
2522b4ed282SShri Abhyankar 
2532b4ed282SShri Abhyankar    Notes:
2542b4ed282SShri Abhyankar    The semismooth jacobian matrix is given by
2552b4ed282SShri Abhyankar    jac = Da + Db*jacfun
2562b4ed282SShri Abhyankar    where Db is the row scaling matrix stored as a vector,
2572b4ed282SShri Abhyankar          Da is the diagonal perturbation matrix stored as a vector
2582b4ed282SShri Abhyankar    and   jacfun is the jacobian of the original nonlinear function.
2592b4ed282SShri Abhyankar */
2602b4ed282SShri Abhyankar #undef __FUNCT__
2611f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian"
2621f79c099SShri Abhyankar PetscErrorCode SNESVIComputeJacobian(SNES snes,Vec X,Mat *jac, Mat *jac_pre, MatStructure *flg,void* jacctx)
2632b4ed282SShri Abhyankar {
2642b4ed282SShri Abhyankar   PetscErrorCode ierr;
2652b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*)snes->data;
2662b4ed282SShri Abhyankar   PetscScalar    *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei;
2672b4ed282SShri Abhyankar   PetscInt       i,nlocal;
2682b4ed282SShri Abhyankar   Vec            F = snes->vec_func;
2692b4ed282SShri Abhyankar 
2702b4ed282SShri Abhyankar   PetscFunctionBegin;
2712b4ed282SShri Abhyankar 
2722b4ed282SShri Abhyankar   ierr = (*vi->computeuserjacobian)(snes,X,jac,jac_pre,flg,jacctx);CHKERRQ(ierr);
2732b4ed282SShri Abhyankar 
2742b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
2752b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
2762b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr);
2772b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr);
2782b4ed282SShri Abhyankar   ierr = VecGetArray(vi->Da,&da);CHKERRQ(ierr);
2792b4ed282SShri Abhyankar   ierr = VecGetArray(vi->Db,&db);CHKERRQ(ierr);
2802b4ed282SShri Abhyankar   ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr);
2812b4ed282SShri Abhyankar 
2822b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2832b4ed282SShri Abhyankar   /* Set the elements of the vector z:
2842b4ed282SShri Abhyankar      z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0)
2852b4ed282SShri Abhyankar      else z[i] = 0
2862b4ed282SShri Abhyankar   */
2872b4ed282SShri Abhyankar   for(i=0;i < nlocal;i++) {
2882b4ed282SShri Abhyankar     da[i] = db[i] = z[i] = 0;
2892b4ed282SShri Abhyankar     if(PetscAbsScalar(f[i]) <= vi->const_tol) {
2902b4ed282SShri Abhyankar       if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) {
2912b4ed282SShri Abhyankar 	da[i] = 1;
2922b4ed282SShri Abhyankar 	z[i]  = 1;
2932b4ed282SShri Abhyankar       }
2942b4ed282SShri Abhyankar       if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) {
2952b4ed282SShri Abhyankar 	db[i] = 1;
2962b4ed282SShri Abhyankar 	z[i]  = 1;
2972b4ed282SShri Abhyankar       }
2982b4ed282SShri Abhyankar     }
2992b4ed282SShri Abhyankar   }
3002b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr);
3012b4ed282SShri Abhyankar   ierr = MatMult(*jac,vi->z,vi->t);CHKERRQ(ierr);
3022b4ed282SShri Abhyankar   ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr);
3032b4ed282SShri Abhyankar   /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */
3042b4ed282SShri Abhyankar   for(i=0;i< nlocal;i++) {
3052b4ed282SShri Abhyankar     /* Free variables */
3062b4ed282SShri Abhyankar     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {
3072b4ed282SShri Abhyankar       da[i] = 0; db[i] = -1;
3082b4ed282SShri Abhyankar     }
3092b4ed282SShri Abhyankar     /* lower bounded variables */
3102b4ed282SShri Abhyankar     else if (u[i] >= PETSC_VI_INF) {
3112b4ed282SShri Abhyankar       if (da[i] >= 1) {
3122b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t[i]);
3132b4ed282SShri Abhyankar 	da[i] = 1/t2 - 1;
3142b4ed282SShri Abhyankar 	db[i] = t[i]/t2 - 1;
3152b4ed282SShri Abhyankar       } else {
3162b4ed282SShri Abhyankar 	t1 = x[i] - l[i];
3172b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(t1,f[i]);
3182b4ed282SShri Abhyankar 	da[i] = t1/t2 - 1;
3192b4ed282SShri Abhyankar 	db[i] = f[i]/t2 - 1;
3202b4ed282SShri Abhyankar       }
3212b4ed282SShri Abhyankar     }
3222b4ed282SShri Abhyankar     /* upper bounded variables */
3232b4ed282SShri Abhyankar     else if (l[i] <= PETSC_VI_NINF) {
3242b4ed282SShri Abhyankar       if (db[i] >= 1) {
3252b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t[i]);
3262b4ed282SShri Abhyankar 	da[i] = -1/t2 -1;
3272b4ed282SShri Abhyankar 	db[i] = -t[i]/t2 - 1;
3282b4ed282SShri Abhyankar       }
3292b4ed282SShri Abhyankar       else {
3302b4ed282SShri Abhyankar 	t1 = u[i] - x[i];
3312b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(t1,f[i]);
3322b4ed282SShri Abhyankar 	da[i] = t1/t2 - 1;
3332b4ed282SShri Abhyankar 	db[i] = -f[i]/t2 - 1;
3342b4ed282SShri Abhyankar       }
3352b4ed282SShri Abhyankar     }
3362b4ed282SShri Abhyankar     /* Fixed variables */
3372b4ed282SShri Abhyankar     else if (l[i] == u[i]) {
3382b4ed282SShri Abhyankar       da[i] = -1;
3392b4ed282SShri Abhyankar       db[i] = 0;
3402b4ed282SShri Abhyankar     }
3412b4ed282SShri Abhyankar     /* Box constrained variables */
3422b4ed282SShri Abhyankar     else {
3432b4ed282SShri Abhyankar       if (db[i] >= 1) {
3442b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t[i]);
3452b4ed282SShri Abhyankar 	ci = 1/t2 + 1;
3462b4ed282SShri Abhyankar 	di = t[i]/t2 + 1;
3472b4ed282SShri Abhyankar       }
3482b4ed282SShri Abhyankar       else {
3492b4ed282SShri Abhyankar 	t1 = x[i] - u[i];
3502b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(t1,f[i]);
3512b4ed282SShri Abhyankar 	ci = t1/t2 + 1;
3522b4ed282SShri Abhyankar 	di = f[i]/t2 + 1;
3532b4ed282SShri Abhyankar       }
3542b4ed282SShri Abhyankar 
3552b4ed282SShri Abhyankar       if (da[i] >= 1) {
3562b4ed282SShri Abhyankar 	t1 = ci + di*t[i];
3572b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t1);
3582b4ed282SShri Abhyankar 	t1 = t1/t2 - 1;
3592b4ed282SShri Abhyankar 	t2 = 1/t2  - 1;
3602b4ed282SShri Abhyankar       }
3612b4ed282SShri Abhyankar       else {
3622b4ed282SShri Abhyankar 	ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr);
3632b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(x[i]-l[i],ei);
3642b4ed282SShri Abhyankar 	t1 = ei/t2 - 1;
3652b4ed282SShri Abhyankar 	t2 = (x[i] - l[i])/t2 - 1;
3662b4ed282SShri Abhyankar       }
3672b4ed282SShri Abhyankar 
3682b4ed282SShri Abhyankar       da[i] = t2 + t1*ci;
3692b4ed282SShri Abhyankar       db[i] = t1*di;
3702b4ed282SShri Abhyankar     }
3712b4ed282SShri Abhyankar   }
3722b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
3732b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
3742b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr);
3752b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr);
3762b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->Da,&da);CHKERRQ(ierr);
3772b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->Db,&db);CHKERRQ(ierr);
3782b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr);
3792b4ed282SShri Abhyankar 
3802b4ed282SShri Abhyankar   /* Do row scaling  and add diagonal perturbation */
3812b4ed282SShri Abhyankar   ierr = MatDiagonalScale(*jac,vi->Db,PETSC_NULL);CHKERRQ(ierr);
3822b4ed282SShri Abhyankar   ierr = MatDiagonalSet(*jac,vi->Da,ADD_VALUES);CHKERRQ(ierr);
3832b4ed282SShri Abhyankar   if (*jac != *jac_pre) { /* If jac and jac_pre are different */
3842b4ed282SShri Abhyankar     ierr = MatDiagonalScale(*jac_pre,vi->Db,PETSC_NULL);
3852b4ed282SShri Abhyankar     ierr = MatDiagonalSet(*jac_pre,vi->Da,ADD_VALUES);CHKERRQ(ierr);
3862b4ed282SShri Abhyankar   }
3872b4ed282SShri Abhyankar 
3882b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3892b4ed282SShri Abhyankar }
3902b4ed282SShri Abhyankar 
3912b4ed282SShri Abhyankar /*
3922b4ed282SShri Abhyankar    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
3932b4ed282SShri Abhyankar 
3942b4ed282SShri Abhyankar    Input Parameters:
3952b4ed282SShri Abhyankar .  phi - semismooth function.
3962b4ed282SShri Abhyankar .  H   - semismooth jacobian
3972b4ed282SShri Abhyankar 
3982b4ed282SShri Abhyankar    Output Parameters:
3992b4ed282SShri Abhyankar .  dpsi - merit function gradient
4002b4ed282SShri Abhyankar 
4012b4ed282SShri Abhyankar    Notes:
4022b4ed282SShri Abhyankar    The merit function gradient is computed as follows
4032b4ed282SShri Abhyankar    dpsi = H^T*phi
4042b4ed282SShri Abhyankar */
4052b4ed282SShri Abhyankar #undef __FUNCT__
4062b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
4072b4ed282SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
4082b4ed282SShri Abhyankar {
4092b4ed282SShri Abhyankar   PetscErrorCode ierr;
4102b4ed282SShri Abhyankar 
4112b4ed282SShri Abhyankar   PetscFunctionBegin;
4122b4ed282SShri Abhyankar   ierr = MatMultTranspose(H,phi,dpsi);
4132b4ed282SShri Abhyankar 
4142b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4152b4ed282SShri Abhyankar }
4162b4ed282SShri Abhyankar 
4172b4ed282SShri Abhyankar /*
4182b4ed282SShri Abhyankar    SNESVISetDescentDirection - Sets the descent direction for the semismooth algorithm
4192b4ed282SShri Abhyankar 
4202b4ed282SShri Abhyankar    Input Parameters:
4212b4ed282SShri Abhyankar .  snes - the SNES context.
4222b4ed282SShri Abhyankar .  dpsi - gradient of the merit function.
4232b4ed282SShri Abhyankar 
4242b4ed282SShri Abhyankar    Output Parameters:
4252b4ed282SShri Abhyankar .  flg  - PETSC_TRUE if the sufficient descent condition is not satisfied.
4262b4ed282SShri Abhyankar 
4272b4ed282SShri Abhyankar    Notes:
4282b4ed282SShri Abhyankar    The condition for the sufficient descent direction is
4292b4ed282SShri Abhyankar         dpsi^T*Y > delta*||Y||^rho
4302b4ed282SShri Abhyankar    where rho, delta are as defined in the SNES_VI structure.
4312b4ed282SShri Abhyankar    If this condition is satisfied then the existing descent direction i.e.
4322b4ed282SShri Abhyankar    the direction given by the linear solve should be used otherwise it should be set to the negative of the merit function gradient i.e -dpsi.
4332b4ed282SShri Abhyankar */
4342b4ed282SShri Abhyankar #undef __FUNCT__
4352b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckDescentDirection"
436ace3abfcSBarry Smith PetscErrorCode SNESVICheckDescentDirection(SNES snes,Vec dpsi, Vec Y,PetscBool* flg)
4372b4ed282SShri Abhyankar {
4382b4ed282SShri Abhyankar   PetscErrorCode  ierr;
4392b4ed282SShri Abhyankar   SNES_VI       *vi = (SNES_VI*)snes->data;
4402b4ed282SShri Abhyankar   PetscScalar     dpsidotY;
4412b4ed282SShri Abhyankar   PetscReal       norm_Y,rhs;
4422b4ed282SShri Abhyankar   const PetscReal rho = vi->rho,delta=vi->delta;
4432b4ed282SShri Abhyankar 
4442b4ed282SShri Abhyankar   PetscFunctionBegin;
4452b4ed282SShri Abhyankar 
4462b4ed282SShri Abhyankar   *flg = PETSC_FALSE;
4472b4ed282SShri Abhyankar   ierr = VecDot(dpsi,Y,&dpsidotY);CHKERRQ(ierr);
4482b4ed282SShri Abhyankar   ierr = VecNormBegin(Y,NORM_2,&norm_Y);CHKERRQ(ierr);
4492b4ed282SShri Abhyankar   ierr = VecNormEnd(Y,NORM_2,&norm_Y);CHKERRQ(ierr);
4502b4ed282SShri Abhyankar 
4512b4ed282SShri Abhyankar   rhs = delta*PetscPowScalar(norm_Y,rho);
4522b4ed282SShri Abhyankar 
4532b4ed282SShri Abhyankar   if (dpsidotY <= rhs) *flg = PETSC_TRUE;
4542b4ed282SShri Abhyankar 
4552b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4562b4ed282SShri Abhyankar }
4572b4ed282SShri Abhyankar 
4582b4ed282SShri Abhyankar /*
4592b4ed282SShri Abhyankar    SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region .
4602b4ed282SShri Abhyankar 
4612b4ed282SShri Abhyankar    Input Parameters:
4622b4ed282SShri Abhyankar .  lb - lower bound.
4632b4ed282SShri Abhyankar .  ub - upper bound.
4642b4ed282SShri Abhyankar 
4652b4ed282SShri Abhyankar    Output Parameters:
4662b4ed282SShri Abhyankar .  X - the readjusted initial guess.
4672b4ed282SShri Abhyankar 
4682b4ed282SShri Abhyankar    Notes:
4692b4ed282SShri Abhyankar    The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i]))
4702b4ed282SShri Abhyankar */
4712b4ed282SShri Abhyankar #undef __FUNCT__
4722b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess"
4732b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub)
4742b4ed282SShri Abhyankar {
4752b4ed282SShri Abhyankar   PetscErrorCode ierr;
4762b4ed282SShri Abhyankar   PetscInt       i,nlocal;
4772b4ed282SShri Abhyankar   PetscScalar    *x,*l,*u;
4782b4ed282SShri Abhyankar 
4792b4ed282SShri Abhyankar   PetscFunctionBegin;
4802b4ed282SShri Abhyankar 
4812b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
4822b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
4832b4ed282SShri Abhyankar   ierr = VecGetArray(lb,&l);CHKERRQ(ierr);
4842b4ed282SShri Abhyankar   ierr = VecGetArray(ub,&u);CHKERRQ(ierr);
4852b4ed282SShri Abhyankar 
4862b4ed282SShri Abhyankar   for(i = 0; i < nlocal; i++) {
4872b4ed282SShri Abhyankar     x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i]));
4882b4ed282SShri Abhyankar   }
4892b4ed282SShri Abhyankar 
4902b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
4912b4ed282SShri Abhyankar   ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr);
4922b4ed282SShri Abhyankar   ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr);
4932b4ed282SShri Abhyankar 
4942b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4952b4ed282SShri Abhyankar }
4962b4ed282SShri Abhyankar 
4972b4ed282SShri Abhyankar /*  --------------------------------------------------------------------
4982b4ed282SShri Abhyankar 
4992b4ed282SShri Abhyankar      This file implements a semismooth truncated Newton method with a line search,
5002b4ed282SShri Abhyankar      for solving a system of nonlinear equations in complementarity form, using the KSP, Vec,
5012b4ed282SShri Abhyankar      and Mat interfaces for linear solvers, vectors, and matrices,
5022b4ed282SShri Abhyankar      respectively.
5032b4ed282SShri Abhyankar 
5042b4ed282SShri Abhyankar      The following basic routines are required for each nonlinear solver:
5052b4ed282SShri Abhyankar           SNESCreate_XXX()          - Creates a nonlinear solver context
5062b4ed282SShri Abhyankar           SNESSetFromOptions_XXX()  - Sets runtime options
5072b4ed282SShri Abhyankar           SNESSolve_XXX()           - Solves the nonlinear system
5082b4ed282SShri Abhyankar           SNESDestroy_XXX()         - Destroys the nonlinear solver context
5092b4ed282SShri Abhyankar      The suffix "_XXX" denotes a particular implementation, in this case
5102b4ed282SShri Abhyankar      we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving
5112b4ed282SShri Abhyankar      systems of nonlinear equations with a line search (LS) method.
5122b4ed282SShri Abhyankar      These routines are actually called via the common user interface
5132b4ed282SShri Abhyankar      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
5142b4ed282SShri Abhyankar      SNESDestroy(), so the application code interface remains identical
5152b4ed282SShri Abhyankar      for all nonlinear solvers.
5162b4ed282SShri Abhyankar 
5172b4ed282SShri Abhyankar      Another key routine is:
5182b4ed282SShri Abhyankar           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
5192b4ed282SShri Abhyankar      by setting data structures and options.   The interface routine SNESSetUp()
5202b4ed282SShri Abhyankar      is not usually called directly by the user, but instead is called by
5212b4ed282SShri Abhyankar      SNESSolve() if necessary.
5222b4ed282SShri Abhyankar 
5232b4ed282SShri Abhyankar      Additional basic routines are:
5242b4ed282SShri Abhyankar           SNESView_XXX()            - Prints details of runtime options that
5252b4ed282SShri Abhyankar                                       have actually been used.
5262b4ed282SShri Abhyankar      These are called by application codes via the interface routines
5272b4ed282SShri Abhyankar      SNESView().
5282b4ed282SShri Abhyankar 
5292b4ed282SShri Abhyankar      The various types of solvers (preconditioners, Krylov subspace methods,
5302b4ed282SShri Abhyankar      nonlinear solvers, timesteppers) are all organized similarly, so the
5312b4ed282SShri Abhyankar      above description applies to these categories also.
5322b4ed282SShri Abhyankar 
5332b4ed282SShri Abhyankar     -------------------------------------------------------------------- */
5342b4ed282SShri Abhyankar /*
535*69c03620SShri Abhyankar    SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton
5362b4ed282SShri Abhyankar    method using a line search.
5372b4ed282SShri Abhyankar 
5382b4ed282SShri Abhyankar    Input Parameters:
5392b4ed282SShri Abhyankar .  snes - the SNES context
5402b4ed282SShri Abhyankar 
5412b4ed282SShri Abhyankar    Output Parameter:
5422b4ed282SShri Abhyankar .  outits - number of iterations until termination
5432b4ed282SShri Abhyankar 
5442b4ed282SShri Abhyankar    Application Interface Routine: SNESSolve()
5452b4ed282SShri Abhyankar 
5462b4ed282SShri Abhyankar    Notes:
5472b4ed282SShri Abhyankar    This implements essentially a semismooth Newton method with a
5482b4ed282SShri Abhyankar    line search.  By default a cubic backtracking line search
5492b4ed282SShri Abhyankar    is employed, as described in the text "Numerical Methods for
5502b4ed282SShri Abhyankar    Unconstrained Optimization and Nonlinear Equations" by Dennis
5512b4ed282SShri Abhyankar    and Schnabel.
5522b4ed282SShri Abhyankar */
5532b4ed282SShri Abhyankar #undef __FUNCT__
554*69c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS"
555*69c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes)
5562b4ed282SShri Abhyankar {
5572b4ed282SShri Abhyankar   SNES_VI          *vi = (SNES_VI*)snes->data;
5582b4ed282SShri Abhyankar   PetscErrorCode     ierr;
5592b4ed282SShri Abhyankar   PetscInt           maxits,i,lits;
560ace3abfcSBarry Smith   PetscBool         lssucceed,changedir;
5612b4ed282SShri Abhyankar   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
5622b4ed282SShri Abhyankar   PetscReal          gnorm,xnorm=0,ynorm;
5632b4ed282SShri Abhyankar   Vec                Y,X,F,G,W;
5642b4ed282SShri Abhyankar   KSPConvergedReason kspreason;
5652b4ed282SShri Abhyankar 
5662b4ed282SShri Abhyankar   PetscFunctionBegin;
5672b4ed282SShri Abhyankar   snes->numFailures            = 0;
5682b4ed282SShri Abhyankar   snes->numLinearSolveFailures = 0;
5692b4ed282SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
5702b4ed282SShri Abhyankar 
5712b4ed282SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
5722b4ed282SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
5732b4ed282SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
5742b4ed282SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
5752b4ed282SShri Abhyankar   G		= snes->work[1];
5762b4ed282SShri Abhyankar   W		= snes->work[2];
5772b4ed282SShri Abhyankar 
5782b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5792b4ed282SShri Abhyankar   snes->iter = 0;
5802b4ed282SShri Abhyankar   snes->norm = 0.0;
5812b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5822b4ed282SShri Abhyankar 
5832b4ed282SShri Abhyankar   ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr);
5842b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
5852b4ed282SShri Abhyankar   if (snes->domainerror) {
5862b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
5872b4ed282SShri Abhyankar     PetscFunctionReturn(0);
5882b4ed282SShri Abhyankar   }
5892b4ed282SShri Abhyankar    /* Compute Merit function */
5902b4ed282SShri Abhyankar   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
5912b4ed282SShri Abhyankar 
5922b4ed282SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
5932b4ed282SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
5942b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5952b4ed282SShri Abhyankar 
5962b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5972b4ed282SShri Abhyankar   snes->norm = vi->phinorm;
5982b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5992b4ed282SShri Abhyankar   SNESLogConvHistory(snes,vi->phinorm,0);
6002b4ed282SShri Abhyankar   SNESMonitor(snes,0,vi->phinorm);
6012b4ed282SShri Abhyankar 
6022b4ed282SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
6032b4ed282SShri Abhyankar   snes->ttol = vi->phinorm*snes->rtol;
6042b4ed282SShri Abhyankar   /* test convergence */
6052b4ed282SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
6062b4ed282SShri Abhyankar   if (snes->reason) PetscFunctionReturn(0);
6072b4ed282SShri Abhyankar 
6082b4ed282SShri Abhyankar   for (i=0; i<maxits; i++) {
6092b4ed282SShri Abhyankar 
6102b4ed282SShri Abhyankar     /* Call general purpose update function */
6112b4ed282SShri Abhyankar     if (snes->ops->update) {
6122b4ed282SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
6132b4ed282SShri Abhyankar     }
6142b4ed282SShri Abhyankar 
6152b4ed282SShri Abhyankar     /* Solve J Y = Phi, where J is the semismooth jacobian */
6162b4ed282SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
6172b4ed282SShri Abhyankar 
6182b4ed282SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
6192b4ed282SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
6202b4ed282SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
6212b4ed282SShri Abhyankar     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
6222b4ed282SShri Abhyankar     ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr);
6232b4ed282SShri Abhyankar     if (kspreason < 0 || changedir) {
6242b4ed282SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
6252b4ed282SShri Abhyankar         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
6262b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
6272b4ed282SShri Abhyankar         break;
6282b4ed282SShri Abhyankar       }
6292b4ed282SShri Abhyankar       ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr);
6302b4ed282SShri Abhyankar     }
6312b4ed282SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
6322b4ed282SShri Abhyankar     snes->linear_its += lits;
6332b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
6342b4ed282SShri Abhyankar     /*
6352b4ed282SShri Abhyankar     if (vi->precheckstep) {
636ace3abfcSBarry Smith       PetscBool changed_y = PETSC_FALSE;
6372b4ed282SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
6382b4ed282SShri Abhyankar     }
6392b4ed282SShri Abhyankar 
6402b4ed282SShri Abhyankar     if (PetscLogPrintInfo){
6412b4ed282SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
6422b4ed282SShri Abhyankar     }
6432b4ed282SShri Abhyankar     */
6442b4ed282SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
6452b4ed282SShri Abhyankar          Y <- X - lambda*Y
6462b4ed282SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
6472b4ed282SShri Abhyankar     */
6482b4ed282SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
6492b4ed282SShri Abhyankar     ynorm = 1; gnorm = vi->phinorm;
6502b4ed282SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
6512b4ed282SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
6522b4ed282SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
6532b4ed282SShri Abhyankar     if (snes->domainerror) {
6542b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
6552b4ed282SShri Abhyankar       PetscFunctionReturn(0);
6562b4ed282SShri Abhyankar     }
6572b4ed282SShri Abhyankar     if (!lssucceed) {
6582b4ed282SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
659ace3abfcSBarry Smith 	PetscBool ismin;
6602b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6612b4ed282SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
6622b4ed282SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
6632b4ed282SShri Abhyankar         break;
6642b4ed282SShri Abhyankar       }
6652b4ed282SShri Abhyankar     }
6662b4ed282SShri Abhyankar     /* Update function and solution vectors */
6672b4ed282SShri Abhyankar     vi->phinorm = gnorm;
6682b4ed282SShri Abhyankar     vi->merit = 0.5*vi->phinorm*vi->phinorm;
6692b4ed282SShri Abhyankar     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
6702b4ed282SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
6712b4ed282SShri Abhyankar     /* Monitor convergence */
6722b4ed282SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
6732b4ed282SShri Abhyankar     snes->iter = i+1;
6742b4ed282SShri Abhyankar     snes->norm = vi->phinorm;
6752b4ed282SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
6762b4ed282SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
6772b4ed282SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
6782b4ed282SShri Abhyankar     /* Test for convergence, xnorm = || X || */
6792b4ed282SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
6802b4ed282SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
6812b4ed282SShri Abhyankar     if (snes->reason) break;
6822b4ed282SShri Abhyankar   }
6832b4ed282SShri Abhyankar   if (i == maxits) {
6842b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
6852b4ed282SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
6862b4ed282SShri Abhyankar   }
6872b4ed282SShri Abhyankar   PetscFunctionReturn(0);
6882b4ed282SShri Abhyankar }
689*69c03620SShri Abhyankar 
690*69c03620SShri Abhyankar #undef __FUNCT__
691*69c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS"
692*69c03620SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscScalar thresh,IS* ISact,IS* ISinact)
693*69c03620SShri Abhyankar {
694*69c03620SShri Abhyankar   PetscErrorCode ierr;
695*69c03620SShri Abhyankar   PetscInt       i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0;
696*69c03620SShri Abhyankar   PetscInt       *idx_act,*idx_inact,i1=0,i2=0;
697*69c03620SShri Abhyankar   PetscScalar    *db;
698*69c03620SShri Abhyankar 
699*69c03620SShri Abhyankar   PetscFunctionBegin;
700*69c03620SShri Abhyankar 
701*69c03620SShri Abhyankar   ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr);
702*69c03620SShri Abhyankar   ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr);
703*69c03620SShri Abhyankar   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
704*69c03620SShri Abhyankar   /* Compute the sizes of the active and inactive sets */
705*69c03620SShri Abhyankar   for (i=0; i < nlocal;i++) {
706*69c03620SShri Abhyankar     if (db[i] <= thresh) nloc_isact++;
707*69c03620SShri Abhyankar     else nloc_isact++;
708*69c03620SShri Abhyankar   }
709*69c03620SShri Abhyankar 
710*69c03620SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
711*69c03620SShri Abhyankar   ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr);
712*69c03620SShri Abhyankar 
713*69c03620SShri Abhyankar   /* Creating the indexing arrays */
714*69c03620SShri Abhyankar   for(i=ilow; i < ihigh; i++) {
715*69c03620SShri Abhyankar     if (db[i] <= thresh) idx_act[i1++] = i;
716*69c03620SShri Abhyankar     else idx_inact[i2++] = i;
717*69c03620SShri Abhyankar   }
718*69c03620SShri Abhyankar 
719*69c03620SShri Abhyankar   /* Create the index sets */
720*69c03620SShri Abhyankar   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr);
721*69c03620SShri Abhyankar   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr);
722*69c03620SShri Abhyankar 
723*69c03620SShri Abhyankar   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
724*69c03620SShri Abhyankar   ierr = PetscFree(idx_act);CHKERRQ(ierr);
725*69c03620SShri Abhyankar   ierr = PetscFree(idx_inact);CHKERRQ(ierr);
726*69c03620SShri Abhyankar   PetscFunctionReturn(0);
727*69c03620SShri Abhyankar }
728*69c03620SShri Abhyankar 
729*69c03620SShri Abhyankar /* Variational Inequality solver using active set method */
730*69c03620SShri Abhyankar #undef __FUNCT__
731*69c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS"
732*69c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes)
733*69c03620SShri Abhyankar {
734*69c03620SShri Abhyankar   SNES_VI          *vi = (SNES_VI*)snes->data;
735*69c03620SShri Abhyankar   PetscErrorCode     ierr;
736*69c03620SShri Abhyankar   PetscInt           maxits,i,lits;
737*69c03620SShri Abhyankar   PetscBool         lssucceed,changedir;
738*69c03620SShri Abhyankar   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
739*69c03620SShri Abhyankar   PetscReal          gnorm,xnorm=0,ynorm;
740*69c03620SShri Abhyankar   Vec                Y,X,F,G,W;
741*69c03620SShri Abhyankar   KSPConvergedReason kspreason;
742*69c03620SShri Abhyankar   IS                 IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
743*69c03620SShri Abhyankar   PetscScalar        thresh,J_norm1;
744*69c03620SShri Abhyankar 
745*69c03620SShri Abhyankar   PetscFunctionBegin;
746*69c03620SShri Abhyankar   snes->numFailures            = 0;
747*69c03620SShri Abhyankar   snes->numLinearSolveFailures = 0;
748*69c03620SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
749*69c03620SShri Abhyankar 
750*69c03620SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
751*69c03620SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
752*69c03620SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
753*69c03620SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
754*69c03620SShri Abhyankar   G		= snes->work[1];
755*69c03620SShri Abhyankar   W		= snes->work[2];
756*69c03620SShri Abhyankar 
757*69c03620SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
758*69c03620SShri Abhyankar   snes->iter = 0;
759*69c03620SShri Abhyankar   snes->norm = 0.0;
760*69c03620SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
761*69c03620SShri Abhyankar 
762*69c03620SShri Abhyankar   ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr);
763*69c03620SShri Abhyankar   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
764*69c03620SShri Abhyankar   if (snes->domainerror) {
765*69c03620SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
766*69c03620SShri Abhyankar     PetscFunctionReturn(0);
767*69c03620SShri Abhyankar   }
768*69c03620SShri Abhyankar    /* Compute Merit function */
769*69c03620SShri Abhyankar   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
770*69c03620SShri Abhyankar 
771*69c03620SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
772*69c03620SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
773*69c03620SShri Abhyankar   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
774*69c03620SShri Abhyankar 
775*69c03620SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
776*69c03620SShri Abhyankar   snes->norm = vi->phinorm;
777*69c03620SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
778*69c03620SShri Abhyankar   SNESLogConvHistory(snes,vi->phinorm,0);
779*69c03620SShri Abhyankar   SNESMonitor(snes,0,vi->phinorm);
780*69c03620SShri Abhyankar 
781*69c03620SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
782*69c03620SShri Abhyankar   snes->ttol = vi->phinorm*snes->rtol;
783*69c03620SShri Abhyankar   /* test convergence */
784*69c03620SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
785*69c03620SShri Abhyankar   if (snes->reason) PetscFunctionReturn(0);
786*69c03620SShri Abhyankar 
787*69c03620SShri Abhyankar   for (i=0; i<maxits; i++) {
788*69c03620SShri Abhyankar 
789*69c03620SShri Abhyankar     /* Call general purpose update function */
790*69c03620SShri Abhyankar     if (snes->ops->update) {
791*69c03620SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
792*69c03620SShri Abhyankar     }
793*69c03620SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
794*69c03620SShri Abhyankar     /* Compute the threshold value for creating active and inactive sets */
795*69c03620SShri Abhyankar     ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr);
796*69c03620SShri Abhyankar     thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1);
797*69c03620SShri Abhyankar     /* Create active and inactive index sets */
798*69c03620SShri Abhyankar     ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr);
799*69c03620SShri Abhyankar     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Active set semismooth algorithm not implemented yet");
800*69c03620SShri Abhyankar     ierr = VecView(vi->Db,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
801*69c03620SShri Abhyankar     ierr = ISView(IS_act,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
802*69c03620SShri Abhyankar     ierr = ISView(IS_inact,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
803*69c03620SShri Abhyankar 
804*69c03620SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
805*69c03620SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
806*69c03620SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
807*69c03620SShri Abhyankar     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
808*69c03620SShri Abhyankar     ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr);
809*69c03620SShri Abhyankar     if (kspreason < 0 || changedir) {
810*69c03620SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
811*69c03620SShri Abhyankar         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
812*69c03620SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
813*69c03620SShri Abhyankar         break;
814*69c03620SShri Abhyankar       }
815*69c03620SShri Abhyankar       ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr);
816*69c03620SShri Abhyankar     }
817*69c03620SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
818*69c03620SShri Abhyankar     snes->linear_its += lits;
819*69c03620SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
820*69c03620SShri Abhyankar     /*
821*69c03620SShri Abhyankar     if (vi->precheckstep) {
822*69c03620SShri Abhyankar       PetscBool changed_y = PETSC_FALSE;
823*69c03620SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
824*69c03620SShri Abhyankar     }
825*69c03620SShri Abhyankar 
826*69c03620SShri Abhyankar     if (PetscLogPrintInfo){
827*69c03620SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
828*69c03620SShri Abhyankar     }
829*69c03620SShri Abhyankar     */
830*69c03620SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
831*69c03620SShri Abhyankar          Y <- X - lambda*Y
832*69c03620SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
833*69c03620SShri Abhyankar     */
834*69c03620SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
835*69c03620SShri Abhyankar     ynorm = 1; gnorm = vi->phinorm;
836*69c03620SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
837*69c03620SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
838*69c03620SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
839*69c03620SShri Abhyankar     if (snes->domainerror) {
840*69c03620SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
841*69c03620SShri Abhyankar       PetscFunctionReturn(0);
842*69c03620SShri Abhyankar     }
843*69c03620SShri Abhyankar     if (!lssucceed) {
844*69c03620SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
845*69c03620SShri Abhyankar 	PetscBool ismin;
846*69c03620SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
847*69c03620SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
848*69c03620SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
849*69c03620SShri Abhyankar         break;
850*69c03620SShri Abhyankar       }
851*69c03620SShri Abhyankar     }
852*69c03620SShri Abhyankar     /* Update function and solution vectors */
853*69c03620SShri Abhyankar     vi->phinorm = gnorm;
854*69c03620SShri Abhyankar     vi->merit = 0.5*vi->phinorm*vi->phinorm;
855*69c03620SShri Abhyankar     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
856*69c03620SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
857*69c03620SShri Abhyankar     /* Monitor convergence */
858*69c03620SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
859*69c03620SShri Abhyankar     snes->iter = i+1;
860*69c03620SShri Abhyankar     snes->norm = vi->phinorm;
861*69c03620SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
862*69c03620SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
863*69c03620SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
864*69c03620SShri Abhyankar     /* Test for convergence, xnorm = || X || */
865*69c03620SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
866*69c03620SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
867*69c03620SShri Abhyankar     if (snes->reason) break;
868*69c03620SShri Abhyankar   }
869*69c03620SShri Abhyankar   if (i == maxits) {
870*69c03620SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
871*69c03620SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
872*69c03620SShri Abhyankar   }
873*69c03620SShri Abhyankar   ierr = ISDestroy(IS_act);CHKERRQ(ierr);
874*69c03620SShri Abhyankar   ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
875*69c03620SShri Abhyankar   PetscFunctionReturn(0);
876*69c03620SShri Abhyankar }
877*69c03620SShri Abhyankar 
8782b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
8792b4ed282SShri Abhyankar /*
8802b4ed282SShri Abhyankar    SNESSetUp_VI - Sets up the internal data structures for the later use
8812b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
8822b4ed282SShri Abhyankar 
8832b4ed282SShri Abhyankar    Input Parameter:
8842b4ed282SShri Abhyankar .  snes - the SNES context
8852b4ed282SShri Abhyankar .  x - the solution vector
8862b4ed282SShri Abhyankar 
8872b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
8882b4ed282SShri Abhyankar 
8892b4ed282SShri Abhyankar    Notes:
8902b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
8912b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
8922b4ed282SShri Abhyankar    the call to SNESSolve().
8932b4ed282SShri Abhyankar  */
8942b4ed282SShri Abhyankar #undef __FUNCT__
8952b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
8962b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
8972b4ed282SShri Abhyankar {
8982b4ed282SShri Abhyankar   PetscErrorCode ierr;
8992b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*) snes->data;
9002b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
9012b4ed282SShri Abhyankar 
9022b4ed282SShri Abhyankar   PetscFunctionBegin;
9032b4ed282SShri Abhyankar   if (!snes->vec_sol_update) {
9042b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
9052b4ed282SShri Abhyankar     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
9062b4ed282SShri Abhyankar   }
9072b4ed282SShri Abhyankar   if (!snes->work) {
9082b4ed282SShri Abhyankar     snes->nwork = 3;
9092b4ed282SShri Abhyankar     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
9102b4ed282SShri Abhyankar     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
9112b4ed282SShri Abhyankar   }
9122b4ed282SShri Abhyankar 
9132b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr);
9142b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr);
9152b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr);
9162b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr);
9172b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
9182b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr);
9192b4ed282SShri Abhyankar 
9202b4ed282SShri Abhyankar   /* If the lower and upper bound on variables are not set, set it to
9212b4ed282SShri Abhyankar      -Inf and Inf */
9222b4ed282SShri Abhyankar   if (!vi->xl && !vi->xu) {
9232b4ed282SShri Abhyankar     vi->usersetxbounds = PETSC_FALSE;
9242b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr);
9252b4ed282SShri Abhyankar     ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr);
9262b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr);
9272b4ed282SShri Abhyankar     ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr);
9282b4ed282SShri Abhyankar   } else {
9292b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
9302b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
9312b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
9322b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
9332b4ed282SShri Abhyankar     if ((i_start[0] != i_start[1]) || (i_start[0] != i_start[2]) || (i_end[0] != i_end[1]) || (i_end[0] != i_end[2]))
9342b4ed282SShri Abhyankar       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Distribution of lower bound, upper bound and the solution vector should be identical across all the processors.");
9352b4ed282SShri Abhyankar   }
9362b4ed282SShri Abhyankar 
9372b4ed282SShri Abhyankar   vi->computeuserfunction = snes->ops->computefunction;
9382b4ed282SShri Abhyankar   vi->computeuserjacobian = snes->ops->computejacobian;
9392b4ed282SShri Abhyankar 
9401f79c099SShri Abhyankar   snes->ops->computefunction = SNESVIComputeFunction;
9411f79c099SShri Abhyankar   snes->ops->computejacobian = SNESVIComputeJacobian;
9422b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9432b4ed282SShri Abhyankar }
9442b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
9452b4ed282SShri Abhyankar /*
9462b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
9472b4ed282SShri Abhyankar    with SNESCreate_VI().
9482b4ed282SShri Abhyankar 
9492b4ed282SShri Abhyankar    Input Parameter:
9502b4ed282SShri Abhyankar .  snes - the SNES context
9512b4ed282SShri Abhyankar 
9522b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
9532b4ed282SShri Abhyankar  */
9542b4ed282SShri Abhyankar #undef __FUNCT__
9552b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
9562b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
9572b4ed282SShri Abhyankar {
9582b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*) snes->data;
9592b4ed282SShri Abhyankar   PetscErrorCode ierr;
9602b4ed282SShri Abhyankar 
9612b4ed282SShri Abhyankar   PetscFunctionBegin;
9622b4ed282SShri Abhyankar   if (snes->vec_sol_update) {
9632b4ed282SShri Abhyankar     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
9642b4ed282SShri Abhyankar     snes->vec_sol_update = PETSC_NULL;
9652b4ed282SShri Abhyankar   }
9662b4ed282SShri Abhyankar   if (snes->nwork) {
9672b4ed282SShri Abhyankar     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
9682b4ed282SShri Abhyankar     snes->nwork = 0;
9692b4ed282SShri Abhyankar     snes->work  = PETSC_NULL;
9702b4ed282SShri Abhyankar   }
9712b4ed282SShri Abhyankar 
9722b4ed282SShri Abhyankar   /* clear vectors */
9732b4ed282SShri Abhyankar   ierr = VecDestroy(vi->phi); CHKERRQ(ierr);
9742b4ed282SShri Abhyankar   ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr);
9752b4ed282SShri Abhyankar   ierr = VecDestroy(vi->Da); CHKERRQ(ierr);
9762b4ed282SShri Abhyankar   ierr = VecDestroy(vi->Db); CHKERRQ(ierr);
9772b4ed282SShri Abhyankar   ierr = VecDestroy(vi->z); CHKERRQ(ierr);
9782b4ed282SShri Abhyankar   ierr = VecDestroy(vi->t); CHKERRQ(ierr);
9792b4ed282SShri Abhyankar   if (!vi->usersetxbounds) {
9802b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xl); CHKERRQ(ierr);
9812b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xu); CHKERRQ(ierr);
9822b4ed282SShri Abhyankar   }
983c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
984c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
9852b4ed282SShri Abhyankar   }
9862b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
9872b4ed282SShri Abhyankar 
9882b4ed282SShri Abhyankar   /* clear composed functions */
9892b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
990c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
9912b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9922b4ed282SShri Abhyankar }
9932b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
9942b4ed282SShri Abhyankar #undef __FUNCT__
9952b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
9962b4ed282SShri Abhyankar 
9972b4ed282SShri Abhyankar /*
9982b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c
9992b4ed282SShri Abhyankar 
10002b4ed282SShri Abhyankar */
1001ace3abfcSBarry Smith PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
10022b4ed282SShri Abhyankar {
10032b4ed282SShri Abhyankar   PetscErrorCode ierr;
10042b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1005ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
10062b4ed282SShri Abhyankar 
10072b4ed282SShri Abhyankar   PetscFunctionBegin;
10082b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
10092b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10102b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
10112b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
10122b4ed282SShri Abhyankar   if (vi->postcheckstep) {
10132b4ed282SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
10142b4ed282SShri Abhyankar   }
10152b4ed282SShri Abhyankar   if (changed_y) {
10162b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
10172b4ed282SShri Abhyankar   }
10182b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
10192b4ed282SShri Abhyankar   if (!snes->domainerror) {
10202b4ed282SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
10212b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
10222b4ed282SShri Abhyankar   }
1023c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1024c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1025c92abb8aSShri Abhyankar   }
10262b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10272b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10282b4ed282SShri Abhyankar }
10292b4ed282SShri Abhyankar 
10302b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10312b4ed282SShri Abhyankar #undef __FUNCT__
10322b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
10332b4ed282SShri Abhyankar 
10342b4ed282SShri Abhyankar /*
10352b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
10362b4ed282SShri Abhyankar */
1037ace3abfcSBarry Smith PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
10382b4ed282SShri Abhyankar {
10392b4ed282SShri Abhyankar   PetscErrorCode ierr;
10402b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1041ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
10422b4ed282SShri Abhyankar 
10432b4ed282SShri Abhyankar   PetscFunctionBegin;
10442b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
10452b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10462b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
10472b4ed282SShri Abhyankar   if (vi->postcheckstep) {
10482b4ed282SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
10492b4ed282SShri Abhyankar   }
10502b4ed282SShri Abhyankar   if (changed_y) {
10512b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
10522b4ed282SShri Abhyankar   }
10532b4ed282SShri Abhyankar 
10542b4ed282SShri Abhyankar   /* don't evaluate function the last time through */
10552b4ed282SShri Abhyankar   if (snes->iter < snes->max_its-1) {
10562b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
10572b4ed282SShri Abhyankar   }
10582b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10592b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10602b4ed282SShri Abhyankar }
10612b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10622b4ed282SShri Abhyankar #undef __FUNCT__
10632b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
10642b4ed282SShri Abhyankar /*
10652b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c
10662b4ed282SShri Abhyankar */
1067ace3abfcSBarry Smith PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
10682b4ed282SShri Abhyankar {
10692b4ed282SShri Abhyankar   /*
10702b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
10712b4ed282SShri Abhyankar      minimization problem:
10722b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
10732b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
10742b4ed282SShri Abhyankar      This function z(x) is same as the merit function for the complementarity problem.
10752b4ed282SShri Abhyankar    */
10762b4ed282SShri Abhyankar 
10772b4ed282SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
10782b4ed282SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
10792b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10802b4ed282SShri Abhyankar   PetscScalar    cinitslope;
10812b4ed282SShri Abhyankar #endif
10822b4ed282SShri Abhyankar   PetscErrorCode ierr;
10832b4ed282SShri Abhyankar   PetscInt       count;
10842b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*)snes->data;
1085ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
10862b4ed282SShri Abhyankar   MPI_Comm       comm;
10872b4ed282SShri Abhyankar 
10882b4ed282SShri Abhyankar   PetscFunctionBegin;
10892b4ed282SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
10902b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10912b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
10922b4ed282SShri Abhyankar 
10932b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
10942b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
1095c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1096c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
10972b4ed282SShri Abhyankar     }
10982b4ed282SShri Abhyankar     *gnorm = fnorm;
10992b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
11002b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
11012b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
11022b4ed282SShri Abhyankar     goto theend1;
11032b4ed282SShri Abhyankar   }
11042b4ed282SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1105c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1106c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
11072b4ed282SShri Abhyankar     }
11082b4ed282SShri Abhyankar     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
11092b4ed282SShri Abhyankar     *ynorm = vi->maxstep;
11102b4ed282SShri Abhyankar   }
11112b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
11122b4ed282SShri Abhyankar   minlambda = vi->minlambda/rellength;
11132b4ed282SShri Abhyankar   /*  ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */
11142b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11152b4ed282SShri Abhyankar   ierr      = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr);
11162b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
11172b4ed282SShri Abhyankar #else
11182b4ed282SShri Abhyankar   ierr      = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr);
11192b4ed282SShri Abhyankar #endif
11202b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
11212b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
11222b4ed282SShri Abhyankar 
11232b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
11242b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
11252b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
11262b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
11272b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
11282b4ed282SShri Abhyankar     goto theend1;
11292b4ed282SShri Abhyankar   }
11302b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
11312b4ed282SShri Abhyankar   if (snes->domainerror) {
11322b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11332b4ed282SShri Abhyankar     PetscFunctionReturn(0);
11342b4ed282SShri Abhyankar   }
11352b4ed282SShri Abhyankar   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
11362b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
11372b4ed282SShri Abhyankar   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
11382b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1139c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1140c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
11412b4ed282SShri Abhyankar     }
11422b4ed282SShri Abhyankar     goto theend1;
11432b4ed282SShri Abhyankar   }
11442b4ed282SShri Abhyankar 
11452b4ed282SShri Abhyankar   /* Fit points with quadratic */
11462b4ed282SShri Abhyankar   lambda     = 1.0;
11472b4ed282SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
11482b4ed282SShri Abhyankar   lambdaprev = lambda;
11492b4ed282SShri Abhyankar   gnormprev  = *gnorm;
11502b4ed282SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
11512b4ed282SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
11522b4ed282SShri Abhyankar   else                         lambda = lambdatemp;
11532b4ed282SShri Abhyankar 
11542b4ed282SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
11552b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
11562b4ed282SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
11572b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
11582b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
11592b4ed282SShri Abhyankar     goto theend1;
11602b4ed282SShri Abhyankar   }
11612b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
11622b4ed282SShri Abhyankar   if (snes->domainerror) {
11632b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11642b4ed282SShri Abhyankar     PetscFunctionReturn(0);
11652b4ed282SShri Abhyankar   }
11662b4ed282SShri Abhyankar   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
11672b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1168c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1169c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
11702b4ed282SShri Abhyankar   }
11712b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1172c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1173c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
11742b4ed282SShri Abhyankar     }
11752b4ed282SShri Abhyankar     goto theend1;
11762b4ed282SShri Abhyankar   }
11772b4ed282SShri Abhyankar 
11782b4ed282SShri Abhyankar   /* Fit points with cubic */
11792b4ed282SShri Abhyankar   count = 1;
11802b4ed282SShri Abhyankar   while (PETSC_TRUE) {
11812b4ed282SShri Abhyankar     if (lambda <= minlambda) {
1182c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1183c92abb8aSShri Abhyankar 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
1184c92abb8aSShri Abhyankar 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr);
11852b4ed282SShri Abhyankar       }
11862b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
11872b4ed282SShri Abhyankar       break;
11882b4ed282SShri Abhyankar     }
11892b4ed282SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
11902b4ed282SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
11912b4ed282SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
11922b4ed282SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
11932b4ed282SShri Abhyankar     d  = b*b - 3*a*initslope;
11942b4ed282SShri Abhyankar     if (d < 0.0) d = 0.0;
11952b4ed282SShri Abhyankar     if (a == 0.0) {
11962b4ed282SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
11972b4ed282SShri Abhyankar     } else {
11982b4ed282SShri Abhyankar       lambdatemp = (-b + sqrt(d))/(3.0*a);
11992b4ed282SShri Abhyankar     }
12002b4ed282SShri Abhyankar     lambdaprev = lambda;
12012b4ed282SShri Abhyankar     gnormprev  = *gnorm;
12022b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
12032b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
12042b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
12052b4ed282SShri Abhyankar     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
12062b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
12072b4ed282SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
12082b4ed282SShri Abhyankar       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
12092b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
12102b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
12112b4ed282SShri Abhyankar       break;
12122b4ed282SShri Abhyankar     }
12132b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12142b4ed282SShri Abhyankar     if (snes->domainerror) {
12152b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
12162b4ed282SShri Abhyankar       PetscFunctionReturn(0);
12172b4ed282SShri Abhyankar     }
12182b4ed282SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12192b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
12202b4ed282SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
1221c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
12222b4ed282SShri Abhyankar 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
12232b4ed282SShri Abhyankar       }
12242b4ed282SShri Abhyankar       break;
12252b4ed282SShri Abhyankar     } else {
1226c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
12272b4ed282SShri Abhyankar         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
12282b4ed282SShri Abhyankar       }
12292b4ed282SShri Abhyankar     }
12302b4ed282SShri Abhyankar     count++;
12312b4ed282SShri Abhyankar   }
12322b4ed282SShri Abhyankar   theend1:
12332b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
12342b4ed282SShri Abhyankar   if (vi->postcheckstep && *flag) {
12352b4ed282SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
12362b4ed282SShri Abhyankar     if (changed_y) {
12372b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
12382b4ed282SShri Abhyankar     }
12392b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
12402b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12412b4ed282SShri Abhyankar       if (snes->domainerror) {
12422b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
12432b4ed282SShri Abhyankar         PetscFunctionReturn(0);
12442b4ed282SShri Abhyankar       }
12452b4ed282SShri Abhyankar       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
12462b4ed282SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
12472b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
12482b4ed282SShri Abhyankar       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
12492b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
12502b4ed282SShri Abhyankar     }
12512b4ed282SShri Abhyankar   }
12522b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
12532b4ed282SShri Abhyankar   PetscFunctionReturn(0);
12542b4ed282SShri Abhyankar }
12552b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
12562b4ed282SShri Abhyankar #undef __FUNCT__
12572b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
12582b4ed282SShri Abhyankar /*
12592b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c
12602b4ed282SShri Abhyankar */
1261ace3abfcSBarry Smith PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
12622b4ed282SShri Abhyankar {
12632b4ed282SShri Abhyankar   /*
12642b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
12652b4ed282SShri Abhyankar      minimization problem:
12662b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
12672b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
12682b4ed282SShri Abhyankar      z(x) is the same as the merit function for the complementarity problem
12692b4ed282SShri Abhyankar    */
12702b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
12712b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
12722b4ed282SShri Abhyankar   PetscScalar    cinitslope;
12732b4ed282SShri Abhyankar #endif
12742b4ed282SShri Abhyankar   PetscErrorCode ierr;
12752b4ed282SShri Abhyankar   PetscInt       count;
12762b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1277ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
12782b4ed282SShri Abhyankar 
12792b4ed282SShri Abhyankar   PetscFunctionBegin;
12802b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
12812b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
12822b4ed282SShri Abhyankar 
12832b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
12842b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
1285c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1286c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
12872b4ed282SShri Abhyankar     }
12882b4ed282SShri Abhyankar     *gnorm = fnorm;
12892b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
12902b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
12912b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
12922b4ed282SShri Abhyankar     goto theend2;
12932b4ed282SShri Abhyankar   }
12942b4ed282SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
12952b4ed282SShri Abhyankar     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
12962b4ed282SShri Abhyankar     *ynorm = vi->maxstep;
12972b4ed282SShri Abhyankar   }
12982b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
12992b4ed282SShri Abhyankar   minlambda = vi->minlambda/rellength;
13002b4ed282SShri Abhyankar   /*  ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */
13012b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
13022b4ed282SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
13032b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
13042b4ed282SShri Abhyankar #else
13052b4ed282SShri Abhyankar   ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr);
13062b4ed282SShri Abhyankar #endif
13072b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
13082b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
13092b4ed282SShri Abhyankar   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
13102b4ed282SShri Abhyankar 
13112b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
13122b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
13132b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
13142b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
13152b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
13162b4ed282SShri Abhyankar     goto theend2;
13172b4ed282SShri Abhyankar   }
13182b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13192b4ed282SShri Abhyankar   if (snes->domainerror) {
13202b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13212b4ed282SShri Abhyankar     PetscFunctionReturn(0);
13222b4ed282SShri Abhyankar   }
13232b4ed282SShri Abhyankar   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13242b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
13252b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1326c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1327c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
13282b4ed282SShri Abhyankar     }
13292b4ed282SShri Abhyankar     goto theend2;
13302b4ed282SShri Abhyankar   }
13312b4ed282SShri Abhyankar 
13322b4ed282SShri Abhyankar   /* Fit points with quadratic */
13332b4ed282SShri Abhyankar   lambda = 1.0;
13342b4ed282SShri Abhyankar   count = 1;
13352b4ed282SShri Abhyankar   while (PETSC_TRUE) {
13362b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
1337c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1338c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
1339c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
13402b4ed282SShri Abhyankar       }
13412b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
13422b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
13432b4ed282SShri Abhyankar       break;
13442b4ed282SShri Abhyankar     }
13452b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
13462b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
13472b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
13482b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
13492b4ed282SShri Abhyankar 
13502b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
13512b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
13522b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
13532b4ed282SShri Abhyankar       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
13542b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
13552b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
13562b4ed282SShri Abhyankar       break;
13572b4ed282SShri Abhyankar     }
13582b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13592b4ed282SShri Abhyankar     if (snes->domainerror) {
13602b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13612b4ed282SShri Abhyankar       PetscFunctionReturn(0);
13622b4ed282SShri Abhyankar     }
13632b4ed282SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13642b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
13652b4ed282SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1366c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1367c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
13682b4ed282SShri Abhyankar       }
13692b4ed282SShri Abhyankar       break;
13702b4ed282SShri Abhyankar     }
13712b4ed282SShri Abhyankar     count++;
13722b4ed282SShri Abhyankar   }
13732b4ed282SShri Abhyankar   theend2:
13742b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
13752b4ed282SShri Abhyankar   if (vi->postcheckstep) {
13762b4ed282SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
13772b4ed282SShri Abhyankar     if (changed_y) {
13782b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
13792b4ed282SShri Abhyankar     }
13802b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
13812b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
13822b4ed282SShri Abhyankar       if (snes->domainerror) {
13832b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13842b4ed282SShri Abhyankar         PetscFunctionReturn(0);
13852b4ed282SShri Abhyankar       }
13862b4ed282SShri Abhyankar       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
13872b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
13882b4ed282SShri Abhyankar       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
13892b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
13902b4ed282SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
13912b4ed282SShri Abhyankar     }
13922b4ed282SShri Abhyankar   }
13932b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13942b4ed282SShri Abhyankar   PetscFunctionReturn(0);
13952b4ed282SShri Abhyankar }
13962b4ed282SShri Abhyankar 
1397ace3abfcSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/
13982b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
13992b4ed282SShri Abhyankar EXTERN_C_BEGIN
14002b4ed282SShri Abhyankar #undef __FUNCT__
14012b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI"
14022b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
14032b4ed282SShri Abhyankar {
14042b4ed282SShri Abhyankar   PetscFunctionBegin;
14052b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->LineSearch = func;
14062b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->lsP        = lsctx;
14072b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14082b4ed282SShri Abhyankar }
14092b4ed282SShri Abhyankar EXTERN_C_END
14102b4ed282SShri Abhyankar 
14112b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
14122b4ed282SShri Abhyankar EXTERN_C_BEGIN
14132b4ed282SShri Abhyankar #undef __FUNCT__
14142b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
1415ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
14162b4ed282SShri Abhyankar {
14172b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
14182b4ed282SShri Abhyankar   PetscErrorCode ierr;
14192b4ed282SShri Abhyankar 
14202b4ed282SShri Abhyankar   PetscFunctionBegin;
1421c92abb8aSShri Abhyankar   if (flg && !vi->lsmonitor) {
1422c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr);
1423c92abb8aSShri Abhyankar   } else if (!flg && vi->lsmonitor) {
1424c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
14252b4ed282SShri Abhyankar   }
14262b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14272b4ed282SShri Abhyankar }
14282b4ed282SShri Abhyankar EXTERN_C_END
14292b4ed282SShri Abhyankar 
14302b4ed282SShri Abhyankar /*
14312b4ed282SShri Abhyankar    SNESView_VI - Prints info from the SNESVI data structure.
14322b4ed282SShri Abhyankar 
14332b4ed282SShri Abhyankar    Input Parameters:
14342b4ed282SShri Abhyankar .  SNES - the SNES context
14352b4ed282SShri Abhyankar .  viewer - visualization context
14362b4ed282SShri Abhyankar 
14372b4ed282SShri Abhyankar    Application Interface Routine: SNESView()
14382b4ed282SShri Abhyankar */
14392b4ed282SShri Abhyankar #undef __FUNCT__
14402b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI"
14412b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
14422b4ed282SShri Abhyankar {
14432b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
14442b4ed282SShri Abhyankar   const char     *cstr;
14452b4ed282SShri Abhyankar   PetscErrorCode ierr;
1446ace3abfcSBarry Smith   PetscBool     iascii;
14472b4ed282SShri Abhyankar 
14482b4ed282SShri Abhyankar   PetscFunctionBegin;
14492b4ed282SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
14502b4ed282SShri Abhyankar   if (iascii) {
14512b4ed282SShri Abhyankar     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
14522b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
14532b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
14542b4ed282SShri Abhyankar     else                                                cstr = "unknown";
14552b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
14562b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
14572b4ed282SShri Abhyankar   } else {
14582b4ed282SShri Abhyankar     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
14592b4ed282SShri Abhyankar   }
14602b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14612b4ed282SShri Abhyankar }
14622b4ed282SShri Abhyankar 
14632b4ed282SShri Abhyankar /*
14642b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
14652b4ed282SShri Abhyankar 
14662b4ed282SShri Abhyankar    Input Parameters:
14672b4ed282SShri Abhyankar .  snes - the SNES context.
14682b4ed282SShri Abhyankar .  xl   - lower bound.
14692b4ed282SShri Abhyankar .  xu   - upper bound.
14702b4ed282SShri Abhyankar 
14712b4ed282SShri Abhyankar    Notes:
14722b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
14732b4ed282SShri Abhyankar    -Infinity and Infinity respectively during SNESSetUp.
14742b4ed282SShri Abhyankar */
14752b4ed282SShri Abhyankar 
14762b4ed282SShri Abhyankar #undef __FUNCT__
14772b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds"
1478*69c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
14792b4ed282SShri Abhyankar {
14802b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
14812b4ed282SShri Abhyankar 
14822b4ed282SShri Abhyankar   PetscFunctionBegin;
14832b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
14842b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
14852b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
14862b4ed282SShri Abhyankar 
14872b4ed282SShri Abhyankar   /* Check if SNESSetFunction is called */
14882b4ed282SShri Abhyankar   if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
14892b4ed282SShri Abhyankar 
14902b4ed282SShri Abhyankar   /* Check if the vector sizes are compatible for lower and upper bounds */
14912b4ed282SShri Abhyankar   if (xl->map->N != snes->vec_func->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xl->map->N,snes->vec_func->map->N);
14922b4ed282SShri Abhyankar   if (xu->map->N != snes->vec_func->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xu->map->N,snes->vec_func->map->N);
14932b4ed282SShri Abhyankar   vi->usersetxbounds = PETSC_TRUE;
14942b4ed282SShri Abhyankar   vi->xl = xl;
14952b4ed282SShri Abhyankar   vi->xu = xu;
14962b4ed282SShri Abhyankar 
14972b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14982b4ed282SShri Abhyankar }
14992b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15002b4ed282SShri Abhyankar /*
15012b4ed282SShri Abhyankar    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
15022b4ed282SShri Abhyankar 
15032b4ed282SShri Abhyankar    Input Parameter:
15042b4ed282SShri Abhyankar .  snes - the SNES context
15052b4ed282SShri Abhyankar 
15062b4ed282SShri Abhyankar    Application Interface Routine: SNESSetFromOptions()
15072b4ed282SShri Abhyankar */
15082b4ed282SShri Abhyankar #undef __FUNCT__
15092b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI"
15102b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
15112b4ed282SShri Abhyankar {
15122b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
15132b4ed282SShri Abhyankar   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1514*69c03620SShri Abhyankar   const char     *vies[] = {"ss","as"};
15152b4ed282SShri Abhyankar   PetscErrorCode ierr;
15162b4ed282SShri Abhyankar   PetscInt       indx;
1517*69c03620SShri Abhyankar   PetscBool     flg,set,flg2;
15182b4ed282SShri Abhyankar 
15192b4ed282SShri Abhyankar   PetscFunctionBegin;
15202b4ed282SShri Abhyankar     ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
15212b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
15222b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
15232b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
15242b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr);
15252b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr);
15262b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
1527acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
15282b4ed282SShri Abhyankar     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
1529*69c03620SShri Abhyankar     ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr);
1530*69c03620SShri Abhyankar     if (flg2) {
1531*69c03620SShri Abhyankar       switch (indx) {
1532*69c03620SShri Abhyankar       case 0:
1533*69c03620SShri Abhyankar 	snes->ops->solve = SNESSolveVI_SS;
1534*69c03620SShri Abhyankar 	break;
1535*69c03620SShri Abhyankar       case 1:
1536*69c03620SShri Abhyankar 	snes->ops->solve = SNESSolveVI_AS;
1537*69c03620SShri Abhyankar 	break;
1538*69c03620SShri Abhyankar       }
1539*69c03620SShri Abhyankar     }
1540c92abb8aSShri Abhyankar     ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
15412b4ed282SShri Abhyankar     if (flg) {
15422b4ed282SShri Abhyankar       switch (indx) {
15432b4ed282SShri Abhyankar       case 0:
15442b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
15452b4ed282SShri Abhyankar         break;
15462b4ed282SShri Abhyankar       case 1:
15472b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
15482b4ed282SShri Abhyankar         break;
15492b4ed282SShri Abhyankar       case 2:
15502b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
15512b4ed282SShri Abhyankar         break;
15522b4ed282SShri Abhyankar       case 3:
15532b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
15542b4ed282SShri Abhyankar         break;
15552b4ed282SShri Abhyankar       }
15562b4ed282SShri Abhyankar     }
15572b4ed282SShri Abhyankar   ierr = PetscOptionsTail();CHKERRQ(ierr);
15582b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15592b4ed282SShri Abhyankar }
15602b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15612b4ed282SShri Abhyankar /*MC
15622b4ed282SShri Abhyankar       SNESVI - Semismooth newton method based nonlinear solver that uses a line search
15632b4ed282SShri Abhyankar 
15642b4ed282SShri Abhyankar    Options Database:
15652b4ed282SShri Abhyankar +   -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search
15662b4ed282SShri Abhyankar .   -snes_vi_alpha <alpha> - Sets alpha
15672b4ed282SShri Abhyankar .   -snes_vi_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
15682b4ed282SShri Abhyankar .   -snes_vi_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
15692b4ed282SShri Abhyankar     -snes_vi_delta <delta> - Sets the fraction used in the descent test.
15702b4ed282SShri Abhyankar     -snes_vi_rho <rho> - Sets the power used in the descent test.
15712b4ed282SShri Abhyankar      For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho
15722b4ed282SShri Abhyankar -   -snes_vi_monitor - print information about progress of line searches
15732b4ed282SShri Abhyankar 
15742b4ed282SShri Abhyankar 
15752b4ed282SShri Abhyankar    Level: beginner
15762b4ed282SShri Abhyankar 
15772b4ed282SShri Abhyankar .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
15782b4ed282SShri Abhyankar            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
15792b4ed282SShri Abhyankar            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
15802b4ed282SShri Abhyankar 
15812b4ed282SShri Abhyankar M*/
15822b4ed282SShri Abhyankar EXTERN_C_BEGIN
15832b4ed282SShri Abhyankar #undef __FUNCT__
15842b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI"
15852b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes)
15862b4ed282SShri Abhyankar {
15872b4ed282SShri Abhyankar   PetscErrorCode ierr;
15882b4ed282SShri Abhyankar   SNES_VI      *vi;
15892b4ed282SShri Abhyankar 
15902b4ed282SShri Abhyankar   PetscFunctionBegin;
15912b4ed282SShri Abhyankar   snes->ops->setup	     = SNESSetUp_VI;
1592*69c03620SShri Abhyankar   snes->ops->solve	     = SNESSolveVI_SS;
15932b4ed282SShri Abhyankar   snes->ops->destroy	     = SNESDestroy_VI;
15942b4ed282SShri Abhyankar   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
15952b4ed282SShri Abhyankar   snes->ops->view            = SNESView_VI;
15962b4ed282SShri Abhyankar   snes->ops->converged       = SNESDefaultConverged_VI;
15972b4ed282SShri Abhyankar 
15982b4ed282SShri Abhyankar   ierr                   = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
15992b4ed282SShri Abhyankar   snes->data    	 = (void*)vi;
16002b4ed282SShri Abhyankar   vi->alpha		 = 1.e-4;
16012b4ed282SShri Abhyankar   vi->maxstep		 = 1.e8;
16022b4ed282SShri Abhyankar   vi->minlambda         = 1.e-12;
16032b4ed282SShri Abhyankar   vi->LineSearch        = SNESLineSearchCubic_VI;
16042b4ed282SShri Abhyankar   vi->lsP               = PETSC_NULL;
16052b4ed282SShri Abhyankar   vi->postcheckstep     = PETSC_NULL;
16062b4ed282SShri Abhyankar   vi->postcheck         = PETSC_NULL;
16072b4ed282SShri Abhyankar   vi->precheckstep      = PETSC_NULL;
16082b4ed282SShri Abhyankar   vi->precheck          = PETSC_NULL;
16092b4ed282SShri Abhyankar   vi->rho               = 2.1;
16102b4ed282SShri Abhyankar   vi->delta             = 1e-10;
16112b4ed282SShri Abhyankar   vi->const_tol         =  2.2204460492503131e-16;
16122b4ed282SShri Abhyankar   vi->computessfunction = ComputeFischerFunction;
16132b4ed282SShri Abhyankar 
16142b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
16152b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
16162b4ed282SShri Abhyankar 
16172b4ed282SShri Abhyankar   PetscFunctionReturn(0);
16182b4ed282SShri Abhyankar }
16192b4ed282SShri Abhyankar EXTERN_C_END
1620