xref: /petsc/src/snes/impls/vi/vi.c (revision 84c105d79988332506723df2f1301d15cbea5100)
12b4ed282SShri Abhyankar #define PETSCSNES_DLL
22b4ed282SShri Abhyankar 
32b4ed282SShri Abhyankar #include "../src/snes/impls/vi/viimpl.h"
4408da460SShri Abhyankar #include "../include/private/kspimpl.h"
52b4ed282SShri Abhyankar 
62b4ed282SShri Abhyankar /*
72b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
82b4ed282SShri 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
92b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
102b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
112b4ed282SShri Abhyankar */
122b4ed282SShri Abhyankar #undef __FUNCT__
132b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
14ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
152b4ed282SShri Abhyankar {
162b4ed282SShri Abhyankar   PetscReal      a1;
172b4ed282SShri Abhyankar   PetscErrorCode ierr;
18ace3abfcSBarry Smith   PetscBool     hastranspose;
192b4ed282SShri Abhyankar 
202b4ed282SShri Abhyankar   PetscFunctionBegin;
212b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
222b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
232b4ed282SShri Abhyankar   if (hastranspose) {
242b4ed282SShri Abhyankar     /* Compute || J^T F|| */
252b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
262b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
272b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
282b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
292b4ed282SShri Abhyankar   } else {
302b4ed282SShri Abhyankar     Vec         work;
312b4ed282SShri Abhyankar     PetscScalar result;
322b4ed282SShri Abhyankar     PetscReal   wnorm;
332b4ed282SShri Abhyankar 
342b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
352b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
362b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
372b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
382b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
392b4ed282SShri Abhyankar     ierr = VecDestroy(work);CHKERRQ(ierr);
402b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
412b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
422b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
432b4ed282SShri Abhyankar   }
442b4ed282SShri Abhyankar   PetscFunctionReturn(0);
452b4ed282SShri Abhyankar }
462b4ed282SShri Abhyankar 
472b4ed282SShri Abhyankar /*
482b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
492b4ed282SShri Abhyankar */
502b4ed282SShri Abhyankar #undef __FUNCT__
512b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
522b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
532b4ed282SShri Abhyankar {
542b4ed282SShri Abhyankar   PetscReal      a1,a2;
552b4ed282SShri Abhyankar   PetscErrorCode ierr;
56ace3abfcSBarry Smith   PetscBool     hastranspose;
572b4ed282SShri Abhyankar 
582b4ed282SShri Abhyankar   PetscFunctionBegin;
592b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
602b4ed282SShri Abhyankar   if (hastranspose) {
612b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
622b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
632b4ed282SShri Abhyankar 
642b4ed282SShri Abhyankar     /* Compute || J^T W|| */
652b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
662b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
672b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
682b4ed282SShri Abhyankar     if (a1 != 0.0) {
692b4ed282SShri Abhyankar       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
702b4ed282SShri Abhyankar     }
712b4ed282SShri Abhyankar   }
722b4ed282SShri Abhyankar   PetscFunctionReturn(0);
732b4ed282SShri Abhyankar }
742b4ed282SShri Abhyankar 
752b4ed282SShri Abhyankar /*
762b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
772b4ed282SShri Abhyankar 
782b4ed282SShri Abhyankar   Notes:
792b4ed282SShri Abhyankar   The convergence criterion currently implemented is
802b4ed282SShri Abhyankar   merit < abstol
812b4ed282SShri Abhyankar   merit < rtol*merit_initial
822b4ed282SShri Abhyankar */
832b4ed282SShri Abhyankar #undef __FUNCT__
842b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
852b4ed282SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy)
862b4ed282SShri Abhyankar {
872b4ed282SShri Abhyankar   PetscErrorCode ierr;
882b4ed282SShri Abhyankar 
892b4ed282SShri Abhyankar   PetscFunctionBegin;
902b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
912b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
922b4ed282SShri Abhyankar 
932b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
942b4ed282SShri Abhyankar 
952b4ed282SShri Abhyankar   if (!it) {
962b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
972b4ed282SShri Abhyankar     snes->ttol = merit*snes->rtol;
982b4ed282SShri Abhyankar   }
992b4ed282SShri Abhyankar   if (merit != merit) {
1002b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
1012b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1022b4ed282SShri Abhyankar   } else if (merit < snes->abstol) {
1032b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr);
1042b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
1052b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
1062b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
1072b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
1082b4ed282SShri Abhyankar   }
1092b4ed282SShri Abhyankar 
1102b4ed282SShri Abhyankar   if (it && !*reason) {
1112b4ed282SShri Abhyankar     if (merit < snes->ttol) {
1122b4ed282SShri Abhyankar       ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr);
1132b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
1142b4ed282SShri Abhyankar     }
1152b4ed282SShri Abhyankar   }
1162b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1172b4ed282SShri Abhyankar }
1182b4ed282SShri Abhyankar 
1192b4ed282SShri Abhyankar /*
1202b4ed282SShri Abhyankar   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
1212b4ed282SShri Abhyankar 
1222b4ed282SShri Abhyankar   Input Parameter:
1232b4ed282SShri Abhyankar . phi - the semismooth function
1242b4ed282SShri Abhyankar 
1252b4ed282SShri Abhyankar   Output Parameter:
1262b4ed282SShri Abhyankar . merit - the merit function
1272b4ed282SShri Abhyankar . phinorm - ||phi||
1282b4ed282SShri Abhyankar 
1292b4ed282SShri Abhyankar   Notes:
1302b4ed282SShri Abhyankar   The merit function for the mixed complementarity problem is defined as
1312b4ed282SShri Abhyankar      merit = 0.5*phi^T*phi
1322b4ed282SShri Abhyankar */
1332b4ed282SShri Abhyankar #undef __FUNCT__
1342b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction"
1352b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
1362b4ed282SShri Abhyankar {
1372b4ed282SShri Abhyankar   PetscErrorCode ierr;
1382b4ed282SShri Abhyankar 
1392b4ed282SShri Abhyankar   PetscFunctionBegin;
1402b4ed282SShri Abhyankar   ierr = VecNormBegin(phi,NORM_2,phinorm);
1412b4ed282SShri Abhyankar   ierr = VecNormEnd(phi,NORM_2,phinorm);
1422b4ed282SShri Abhyankar 
1432b4ed282SShri Abhyankar   *merit = 0.5*(*phinorm)*(*phinorm);
1442b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1452b4ed282SShri Abhyankar }
1462b4ed282SShri Abhyankar 
1472b4ed282SShri Abhyankar /*
1482b4ed282SShri Abhyankar   ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation.
1492b4ed282SShri Abhyankar 
1502b4ed282SShri Abhyankar   Notes:
1512b4ed282SShri Abhyankar   The Fischer-Burmeister function is defined as
1522b4ed282SShri Abhyankar        ff(a,b) := sqrt(a*a + b*b) - a - b
1532b4ed282SShri Abhyankar   and is used reformulate a complementarity equation  as a semismooth equation.
1542b4ed282SShri Abhyankar */
1552b4ed282SShri Abhyankar 
1562b4ed282SShri Abhyankar #undef __FUNCT__
1572b4ed282SShri Abhyankar #define __FUNCT__ "ComputeFischerFunction"
1582b4ed282SShri Abhyankar static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff)
1592b4ed282SShri Abhyankar {
1602b4ed282SShri Abhyankar   PetscFunctionBegin;
1612b4ed282SShri Abhyankar   *ff = sqrt(a*a + b*b) - a - b;
1622b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1632b4ed282SShri Abhyankar }
1642b4ed282SShri Abhyankar 
1652b4ed282SShri Abhyankar /*
1661f79c099SShri Abhyankar    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
1672b4ed282SShri Abhyankar 
1682b4ed282SShri Abhyankar    Input Parameters:
1692b4ed282SShri Abhyankar .  snes - the SNES context
1702b4ed282SShri Abhyankar .  x - current iterate
1712b4ed282SShri Abhyankar .  functx - user defined function context
1722b4ed282SShri Abhyankar 
1732b4ed282SShri Abhyankar    Output Parameters:
1742b4ed282SShri Abhyankar .  phi - Semismooth function
1752b4ed282SShri Abhyankar 
1762b4ed282SShri Abhyankar    Notes:
1772b4ed282SShri Abhyankar    The result of this function is done by cases:
1782b4ed282SShri Abhyankar +  l[i] == -infinity, u[i] == infinity  -- phi[i] = -f[i]
1792b4ed282SShri Abhyankar .  l[i] == -infinity, u[i] finite       -- phi[i] = ss(u[i]-x[i], -f[i])
1802b4ed282SShri Abhyankar .  l[i] finite,       u[i] == infinity  -- phi[i] = ss(x[i]-l[i],  f[i])
1812b4ed282SShri Abhyankar .  l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u]))
1822b4ed282SShri Abhyankar -  otherwise l[i] == u[i] -- phi[i] = l[i] - x[i]
1832b4ed282SShri Abhyankar    ss is the semismoothing function used to reformulate the nonlinear equation in complementarity
1842b4ed282SShri Abhyankar    form to semismooth form
1852b4ed282SShri Abhyankar 
1862b4ed282SShri Abhyankar */
1872b4ed282SShri Abhyankar #undef __FUNCT__
1881f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction"
1891f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
1902b4ed282SShri Abhyankar {
1912b4ed282SShri Abhyankar   PetscErrorCode  ierr;
1922b4ed282SShri Abhyankar   SNES_VI       *vi = (SNES_VI*)snes->data;
1932b4ed282SShri Abhyankar   Vec             Xl = vi->xl,Xu = vi->xu,F = snes->vec_func;
1942b4ed282SShri Abhyankar   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u,t;
1952b4ed282SShri Abhyankar   PetscInt        i,nlocal;
1962b4ed282SShri Abhyankar 
1972b4ed282SShri Abhyankar   PetscFunctionBegin;
1982b4ed282SShri Abhyankar   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
1992b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2002b4ed282SShri Abhyankar 
2012b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
2022b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
2032b4ed282SShri Abhyankar   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
2042b4ed282SShri Abhyankar   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
2052b4ed282SShri Abhyankar   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
2062b4ed282SShri Abhyankar 
2072b4ed282SShri Abhyankar   for (i=0;i < nlocal;i++) {
2082b4ed282SShri Abhyankar     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {
2092b4ed282SShri Abhyankar       phi_arr[i] = -f_arr[i];
2102b4ed282SShri Abhyankar     }
2112b4ed282SShri Abhyankar     else if (l[i] <= PETSC_VI_NINF) {
2122b4ed282SShri Abhyankar       t = u[i] - x_arr[i];
2132b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr);
2142b4ed282SShri Abhyankar       phi_arr[i] = -phi_arr[i];
2152b4ed282SShri Abhyankar     }
2162b4ed282SShri Abhyankar     else if (u[i] >= PETSC_VI_INF) {
2172b4ed282SShri Abhyankar       t = x_arr[i] - l[i];
2182b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr);
2192b4ed282SShri Abhyankar     }
2202b4ed282SShri Abhyankar     else if (l[i] == u[i]) {
2212b4ed282SShri Abhyankar       phi_arr[i] = l[i] - x_arr[i];
2222b4ed282SShri Abhyankar     }
2232b4ed282SShri Abhyankar     else {
2242b4ed282SShri Abhyankar       t = u[i] - x_arr[i];
2252b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);
2262b4ed282SShri Abhyankar       t = x_arr[i] - l[i];
2272b4ed282SShri Abhyankar       ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]);
2282b4ed282SShri Abhyankar     }
2292b4ed282SShri Abhyankar   }
2302b4ed282SShri Abhyankar 
2312b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
2322b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
2332b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
2342b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
2352b4ed282SShri Abhyankar   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
2362b4ed282SShri Abhyankar 
2372b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2382b4ed282SShri Abhyankar }
2392b4ed282SShri Abhyankar 
2402b4ed282SShri Abhyankar /*
241a79edbefSShri Abhyankar    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
242a79edbefSShri Abhyankar                                           the semismooth jacobian.
2432b4ed282SShri Abhyankar */
2442b4ed282SShri Abhyankar #undef __FUNCT__
245a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
246a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
2472b4ed282SShri Abhyankar {
2482b4ed282SShri Abhyankar   PetscErrorCode ierr;
2492b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*)snes->data;
2502b4ed282SShri Abhyankar   PetscScalar    *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei;
2512b4ed282SShri Abhyankar   PetscInt       i,nlocal;
2522b4ed282SShri Abhyankar 
2532b4ed282SShri Abhyankar   PetscFunctionBegin;
2542b4ed282SShri Abhyankar 
2552b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
2562b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
2572b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr);
2582b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr);
259a79edbefSShri Abhyankar   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
260a79edbefSShri Abhyankar   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
2612b4ed282SShri Abhyankar   ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr);
2622b4ed282SShri Abhyankar 
2632b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2642b4ed282SShri Abhyankar   /* Set the elements of the vector z:
2652b4ed282SShri Abhyankar      z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0)
2662b4ed282SShri Abhyankar      else z[i] = 0
2672b4ed282SShri Abhyankar   */
2682b4ed282SShri Abhyankar   for(i=0;i < nlocal;i++) {
2692b4ed282SShri Abhyankar     da[i] = db[i] = z[i] = 0;
2702b4ed282SShri Abhyankar     if(PetscAbsScalar(f[i]) <= vi->const_tol) {
2712b4ed282SShri Abhyankar       if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) {
2722b4ed282SShri Abhyankar 	da[i] = 1;
2732b4ed282SShri Abhyankar 	z[i]  = 1;
2742b4ed282SShri Abhyankar       }
2752b4ed282SShri Abhyankar       if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) {
2762b4ed282SShri Abhyankar 	db[i] = 1;
2772b4ed282SShri Abhyankar 	z[i]  = 1;
2782b4ed282SShri Abhyankar       }
2792b4ed282SShri Abhyankar     }
2802b4ed282SShri Abhyankar   }
2812b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr);
282a79edbefSShri Abhyankar   ierr = MatMult(jac,vi->z,vi->t);CHKERRQ(ierr);
2832b4ed282SShri Abhyankar   ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr);
2842b4ed282SShri Abhyankar   /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */
2852b4ed282SShri Abhyankar   for(i=0;i< nlocal;i++) {
2862b4ed282SShri Abhyankar     /* Free variables */
2872b4ed282SShri Abhyankar     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {
2882b4ed282SShri Abhyankar       da[i] = 0; db[i] = -1;
2892b4ed282SShri Abhyankar     }
2902b4ed282SShri Abhyankar     /* lower bounded variables */
2912b4ed282SShri Abhyankar     else if (u[i] >= PETSC_VI_INF) {
2922b4ed282SShri Abhyankar       if (da[i] >= 1) {
2932b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t[i]);
2942b4ed282SShri Abhyankar 	da[i] = 1/t2 - 1;
2952b4ed282SShri Abhyankar 	db[i] = t[i]/t2 - 1;
2962b4ed282SShri Abhyankar       } else {
2972b4ed282SShri Abhyankar 	t1 = x[i] - l[i];
2982b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(t1,f[i]);
2992b4ed282SShri Abhyankar 	da[i] = t1/t2 - 1;
3002b4ed282SShri Abhyankar 	db[i] = f[i]/t2 - 1;
3012b4ed282SShri Abhyankar       }
3022b4ed282SShri Abhyankar     }
3032b4ed282SShri Abhyankar     /* upper bounded variables */
3042b4ed282SShri Abhyankar     else if (l[i] <= PETSC_VI_NINF) {
3052b4ed282SShri Abhyankar       if (db[i] >= 1) {
3062b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t[i]);
3072b4ed282SShri Abhyankar 	da[i] = -1/t2 -1;
3082b4ed282SShri Abhyankar 	db[i] = -t[i]/t2 - 1;
3092b4ed282SShri Abhyankar       }
3102b4ed282SShri Abhyankar       else {
3112b4ed282SShri Abhyankar 	t1 = u[i] - x[i];
3122b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(t1,f[i]);
3132b4ed282SShri Abhyankar 	da[i] = t1/t2 - 1;
3142b4ed282SShri Abhyankar 	db[i] = -f[i]/t2 - 1;
3152b4ed282SShri Abhyankar       }
3162b4ed282SShri Abhyankar     }
3172b4ed282SShri Abhyankar     /* Fixed variables */
3182b4ed282SShri Abhyankar     else if (l[i] == u[i]) {
3192b4ed282SShri Abhyankar       da[i] = -1;
3202b4ed282SShri Abhyankar       db[i] = 0;
3212b4ed282SShri Abhyankar     }
3222b4ed282SShri Abhyankar     /* Box constrained variables */
3232b4ed282SShri Abhyankar     else {
3242b4ed282SShri Abhyankar       if (db[i] >= 1) {
3252b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t[i]);
3262b4ed282SShri Abhyankar 	ci = 1/t2 + 1;
3272b4ed282SShri Abhyankar 	di = t[i]/t2 + 1;
3282b4ed282SShri Abhyankar       }
3292b4ed282SShri Abhyankar       else {
3302b4ed282SShri Abhyankar 	t1 = x[i] - u[i];
3312b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(t1,f[i]);
3322b4ed282SShri Abhyankar 	ci = t1/t2 + 1;
3332b4ed282SShri Abhyankar 	di = f[i]/t2 + 1;
3342b4ed282SShri Abhyankar       }
3352b4ed282SShri Abhyankar 
3362b4ed282SShri Abhyankar       if (da[i] >= 1) {
3372b4ed282SShri Abhyankar 	t1 = ci + di*t[i];
3382b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(1,t1);
3392b4ed282SShri Abhyankar 	t1 = t1/t2 - 1;
3402b4ed282SShri Abhyankar 	t2 = 1/t2  - 1;
3412b4ed282SShri Abhyankar       }
3422b4ed282SShri Abhyankar       else {
3432b4ed282SShri Abhyankar 	ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr);
3442b4ed282SShri Abhyankar 	t2 = PetscScalarNorm(x[i]-l[i],ei);
3452b4ed282SShri Abhyankar 	t1 = ei/t2 - 1;
3462b4ed282SShri Abhyankar 	t2 = (x[i] - l[i])/t2 - 1;
3472b4ed282SShri Abhyankar       }
3482b4ed282SShri Abhyankar 
3492b4ed282SShri Abhyankar       da[i] = t2 + t1*ci;
3502b4ed282SShri Abhyankar       db[i] = t1*di;
3512b4ed282SShri Abhyankar     }
3522b4ed282SShri Abhyankar   }
3532b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
3542b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
3552b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr);
3562b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr);
357a79edbefSShri Abhyankar   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
358a79edbefSShri Abhyankar   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
3592b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr);
3602b4ed282SShri Abhyankar 
361a79edbefSShri Abhyankar   PetscFunctionReturn(0);
362a79edbefSShri Abhyankar }
363a79edbefSShri Abhyankar 
364a79edbefSShri Abhyankar /*
365a79edbefSShri 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.
366a79edbefSShri Abhyankar 
367a79edbefSShri Abhyankar    Input Parameters:
368a79edbefSShri Abhyankar .  Da       - Diagonal shift vector for the semismooth jacobian.
369a79edbefSShri Abhyankar .  Db       - Row scaling vector for the semismooth jacobian.
370a79edbefSShri Abhyankar 
371a79edbefSShri Abhyankar    Output Parameters:
372a79edbefSShri Abhyankar .  jac      - semismooth jacobian
373a79edbefSShri Abhyankar .  jac_pre  - optional preconditioning matrix
374a79edbefSShri Abhyankar 
375a79edbefSShri Abhyankar    Notes:
376a79edbefSShri Abhyankar    The semismooth jacobian matrix is given by
377a79edbefSShri Abhyankar    jac = Da + Db*jacfun
378a79edbefSShri Abhyankar    where Db is the row scaling matrix stored as a vector,
379a79edbefSShri Abhyankar          Da is the diagonal perturbation matrix stored as a vector
380a79edbefSShri Abhyankar    and   jacfun is the jacobian of the original nonlinear function.
381a79edbefSShri Abhyankar */
382a79edbefSShri Abhyankar #undef __FUNCT__
383a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian"
384a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
385a79edbefSShri Abhyankar {
386a79edbefSShri Abhyankar   PetscErrorCode ierr;
387a79edbefSShri Abhyankar 
3882b4ed282SShri Abhyankar   /* Do row scaling  and add diagonal perturbation */
389a79edbefSShri Abhyankar   ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr);
390a79edbefSShri Abhyankar   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
391a79edbefSShri Abhyankar   if (jac != jac_pre) { /* If jac and jac_pre are different */
392a79edbefSShri Abhyankar     ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL);
393a79edbefSShri Abhyankar     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
3942b4ed282SShri Abhyankar   }
3952b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3962b4ed282SShri Abhyankar }
3972b4ed282SShri Abhyankar 
3982b4ed282SShri Abhyankar /*
3992b4ed282SShri Abhyankar    SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region .
4002b4ed282SShri Abhyankar 
4012b4ed282SShri Abhyankar    Input Parameters:
4022b4ed282SShri Abhyankar .  lb - lower bound.
4032b4ed282SShri Abhyankar .  ub - upper bound.
4042b4ed282SShri Abhyankar 
4052b4ed282SShri Abhyankar    Output Parameters:
4062b4ed282SShri Abhyankar .  X - the readjusted initial guess.
4072b4ed282SShri Abhyankar 
4082b4ed282SShri Abhyankar    Notes:
4092b4ed282SShri Abhyankar    The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i]))
4102b4ed282SShri Abhyankar */
4112b4ed282SShri Abhyankar #undef __FUNCT__
4122b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess"
4132b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub)
4142b4ed282SShri Abhyankar {
4152b4ed282SShri Abhyankar   PetscErrorCode ierr;
4162b4ed282SShri Abhyankar   PetscInt       i,nlocal;
4172b4ed282SShri Abhyankar   PetscScalar    *x,*l,*u;
4182b4ed282SShri Abhyankar 
4192b4ed282SShri Abhyankar   PetscFunctionBegin;
4202b4ed282SShri Abhyankar 
4212b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
4222b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
4232b4ed282SShri Abhyankar   ierr = VecGetArray(lb,&l);CHKERRQ(ierr);
4242b4ed282SShri Abhyankar   ierr = VecGetArray(ub,&u);CHKERRQ(ierr);
4252b4ed282SShri Abhyankar 
4262b4ed282SShri Abhyankar   for(i = 0; i < nlocal; i++) {
4272b4ed282SShri Abhyankar     x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i]));
4282b4ed282SShri Abhyankar   }
4292b4ed282SShri Abhyankar 
4302b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
4312b4ed282SShri Abhyankar   ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr);
4322b4ed282SShri Abhyankar   ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr);
4332b4ed282SShri Abhyankar 
4342b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4352b4ed282SShri Abhyankar }
4362b4ed282SShri Abhyankar 
4372b4ed282SShri Abhyankar /*  --------------------------------------------------------------------
4382b4ed282SShri Abhyankar 
4392b4ed282SShri Abhyankar      This file implements a semismooth truncated Newton method with a line search,
4402b4ed282SShri Abhyankar      for solving a system of nonlinear equations in complementarity form, using the KSP, Vec,
4412b4ed282SShri Abhyankar      and Mat interfaces for linear solvers, vectors, and matrices,
4422b4ed282SShri Abhyankar      respectively.
4432b4ed282SShri Abhyankar 
4442b4ed282SShri Abhyankar      The following basic routines are required for each nonlinear solver:
4452b4ed282SShri Abhyankar           SNESCreate_XXX()          - Creates a nonlinear solver context
4462b4ed282SShri Abhyankar           SNESSetFromOptions_XXX()  - Sets runtime options
4472b4ed282SShri Abhyankar           SNESSolve_XXX()           - Solves the nonlinear system
4482b4ed282SShri Abhyankar           SNESDestroy_XXX()         - Destroys the nonlinear solver context
4492b4ed282SShri Abhyankar      The suffix "_XXX" denotes a particular implementation, in this case
4502b4ed282SShri Abhyankar      we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving
4512b4ed282SShri Abhyankar      systems of nonlinear equations with a line search (LS) method.
4522b4ed282SShri Abhyankar      These routines are actually called via the common user interface
4532b4ed282SShri Abhyankar      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
4542b4ed282SShri Abhyankar      SNESDestroy(), so the application code interface remains identical
4552b4ed282SShri Abhyankar      for all nonlinear solvers.
4562b4ed282SShri Abhyankar 
4572b4ed282SShri Abhyankar      Another key routine is:
4582b4ed282SShri Abhyankar           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
4592b4ed282SShri Abhyankar      by setting data structures and options.   The interface routine SNESSetUp()
4602b4ed282SShri Abhyankar      is not usually called directly by the user, but instead is called by
4612b4ed282SShri Abhyankar      SNESSolve() if necessary.
4622b4ed282SShri Abhyankar 
4632b4ed282SShri Abhyankar      Additional basic routines are:
4642b4ed282SShri Abhyankar           SNESView_XXX()            - Prints details of runtime options that
4652b4ed282SShri Abhyankar                                       have actually been used.
4662b4ed282SShri Abhyankar      These are called by application codes via the interface routines
4672b4ed282SShri Abhyankar      SNESView().
4682b4ed282SShri Abhyankar 
4692b4ed282SShri Abhyankar      The various types of solvers (preconditioners, Krylov subspace methods,
4702b4ed282SShri Abhyankar      nonlinear solvers, timesteppers) are all organized similarly, so the
4712b4ed282SShri Abhyankar      above description applies to these categories also.
4722b4ed282SShri Abhyankar 
4732b4ed282SShri Abhyankar     -------------------------------------------------------------------- */
4742b4ed282SShri Abhyankar /*
47569c03620SShri Abhyankar    SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton
4762b4ed282SShri Abhyankar    method using a line search.
4772b4ed282SShri Abhyankar 
4782b4ed282SShri Abhyankar    Input Parameters:
4792b4ed282SShri Abhyankar .  snes - the SNES context
4802b4ed282SShri Abhyankar 
4812b4ed282SShri Abhyankar    Output Parameter:
4822b4ed282SShri Abhyankar .  outits - number of iterations until termination
4832b4ed282SShri Abhyankar 
4842b4ed282SShri Abhyankar    Application Interface Routine: SNESSolve()
4852b4ed282SShri Abhyankar 
4862b4ed282SShri Abhyankar    Notes:
4872b4ed282SShri Abhyankar    This implements essentially a semismooth Newton method with a
4882b4ed282SShri Abhyankar    line search.  By default a cubic backtracking line search
4892b4ed282SShri Abhyankar    is employed, as described in the text "Numerical Methods for
4902b4ed282SShri Abhyankar    Unconstrained Optimization and Nonlinear Equations" by Dennis
4912b4ed282SShri Abhyankar    and Schnabel.
4922b4ed282SShri Abhyankar */
4932b4ed282SShri Abhyankar #undef __FUNCT__
49469c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS"
49569c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes)
4962b4ed282SShri Abhyankar {
4972b4ed282SShri Abhyankar   SNES_VI            *vi = (SNES_VI*)snes->data;
4982b4ed282SShri Abhyankar   PetscErrorCode     ierr;
4992b4ed282SShri Abhyankar   PetscInt           maxits,i,lits;
5003b336b49SShri Abhyankar   PetscBool          lssucceed;
5012b4ed282SShri Abhyankar   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
5022b4ed282SShri Abhyankar   PetscReal          gnorm,xnorm=0,ynorm;
5032b4ed282SShri Abhyankar   Vec                Y,X,F,G,W;
5042b4ed282SShri Abhyankar   KSPConvergedReason kspreason;
5052b4ed282SShri Abhyankar 
5062b4ed282SShri Abhyankar   PetscFunctionBegin;
5072b4ed282SShri Abhyankar   snes->numFailures            = 0;
5082b4ed282SShri Abhyankar   snes->numLinearSolveFailures = 0;
5092b4ed282SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
5102b4ed282SShri Abhyankar 
5112b4ed282SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
5122b4ed282SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
5132b4ed282SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
5142b4ed282SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
5152b4ed282SShri Abhyankar   G		= snes->work[1];
5162b4ed282SShri Abhyankar   W		= snes->work[2];
5172b4ed282SShri Abhyankar 
5182b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5192b4ed282SShri Abhyankar   snes->iter = 0;
5202b4ed282SShri Abhyankar   snes->norm = 0.0;
5212b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5222b4ed282SShri Abhyankar 
5232b4ed282SShri Abhyankar   ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr);
5242b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
5252b4ed282SShri Abhyankar   if (snes->domainerror) {
5262b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
5272b4ed282SShri Abhyankar     PetscFunctionReturn(0);
5282b4ed282SShri Abhyankar   }
5292b4ed282SShri Abhyankar    /* Compute Merit function */
5302b4ed282SShri Abhyankar   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
5312b4ed282SShri Abhyankar 
5322b4ed282SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
5332b4ed282SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
5342b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5352b4ed282SShri Abhyankar 
5362b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5372b4ed282SShri Abhyankar   snes->norm = vi->phinorm;
5382b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5392b4ed282SShri Abhyankar   SNESLogConvHistory(snes,vi->phinorm,0);
5402b4ed282SShri Abhyankar   SNESMonitor(snes,0,vi->phinorm);
5412b4ed282SShri Abhyankar 
5422b4ed282SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
5432b4ed282SShri Abhyankar   snes->ttol = vi->phinorm*snes->rtol;
5442b4ed282SShri Abhyankar   /* test convergence */
5452b4ed282SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
5462b4ed282SShri Abhyankar   if (snes->reason) PetscFunctionReturn(0);
5472b4ed282SShri Abhyankar 
5482b4ed282SShri Abhyankar   for (i=0; i<maxits; i++) {
5492b4ed282SShri Abhyankar 
5502b4ed282SShri Abhyankar     /* Call general purpose update function */
5512b4ed282SShri Abhyankar     if (snes->ops->update) {
5522b4ed282SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
5532b4ed282SShri Abhyankar     }
5542b4ed282SShri Abhyankar 
5552b4ed282SShri Abhyankar     /* Solve J Y = Phi, where J is the semismooth jacobian */
556a79edbefSShri Abhyankar     /* Get the nonlinear function jacobian */
5572b4ed282SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
558a79edbefSShri Abhyankar     /* Get the diagonal shift and row scaling vectors */
559a79edbefSShri Abhyankar     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
560a79edbefSShri Abhyankar     /* Compute the semismooth jacobian */
561a79edbefSShri Abhyankar     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
5622b4ed282SShri Abhyankar 
5632b4ed282SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
5642b4ed282SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
5652b4ed282SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
5663b336b49SShri Abhyankar 
5673b336b49SShri Abhyankar     if (kspreason < 0) {
5682b4ed282SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
5692b4ed282SShri 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);
5702b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
5712b4ed282SShri Abhyankar         break;
5722b4ed282SShri Abhyankar       }
5732b4ed282SShri Abhyankar     }
5742b4ed282SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
5752b4ed282SShri Abhyankar     snes->linear_its += lits;
5762b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
5772b4ed282SShri Abhyankar     /*
5782b4ed282SShri Abhyankar     if (vi->precheckstep) {
579ace3abfcSBarry Smith       PetscBool changed_y = PETSC_FALSE;
5802b4ed282SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
5812b4ed282SShri Abhyankar     }
5822b4ed282SShri Abhyankar 
5832b4ed282SShri Abhyankar     if (PetscLogPrintInfo){
5842b4ed282SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
5852b4ed282SShri Abhyankar     }
5862b4ed282SShri Abhyankar     */
5872b4ed282SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
5882b4ed282SShri Abhyankar          Y <- X - lambda*Y
5892b4ed282SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
5902b4ed282SShri Abhyankar     */
5912b4ed282SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
5922b4ed282SShri Abhyankar     ynorm = 1; gnorm = vi->phinorm;
5932b4ed282SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
5942b4ed282SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
5952b4ed282SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
5962b4ed282SShri Abhyankar     if (snes->domainerror) {
5972b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
5982b4ed282SShri Abhyankar       PetscFunctionReturn(0);
5992b4ed282SShri Abhyankar     }
6002b4ed282SShri Abhyankar     if (!lssucceed) {
6012b4ed282SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
602ace3abfcSBarry Smith 	PetscBool ismin;
6032b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6042b4ed282SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
6052b4ed282SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
6062b4ed282SShri Abhyankar         break;
6072b4ed282SShri Abhyankar       }
6082b4ed282SShri Abhyankar     }
6092b4ed282SShri Abhyankar     /* Update function and solution vectors */
6102b4ed282SShri Abhyankar     vi->phinorm = gnorm;
6112b4ed282SShri Abhyankar     vi->merit = 0.5*vi->phinorm*vi->phinorm;
6122b4ed282SShri Abhyankar     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
6132b4ed282SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
6142b4ed282SShri Abhyankar     /* Monitor convergence */
6152b4ed282SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
6162b4ed282SShri Abhyankar     snes->iter = i+1;
6172b4ed282SShri Abhyankar     snes->norm = vi->phinorm;
6182b4ed282SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
6192b4ed282SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
6202b4ed282SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
6212b4ed282SShri Abhyankar     /* Test for convergence, xnorm = || X || */
6222b4ed282SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
6232b4ed282SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
6242b4ed282SShri Abhyankar     if (snes->reason) break;
6252b4ed282SShri Abhyankar   }
6262b4ed282SShri Abhyankar   if (i == maxits) {
6272b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
6282b4ed282SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
6292b4ed282SShri Abhyankar   }
6302b4ed282SShri Abhyankar   PetscFunctionReturn(0);
6312b4ed282SShri Abhyankar }
63269c03620SShri Abhyankar 
63369c03620SShri Abhyankar #undef __FUNCT__
63469c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS"
635a79edbefSShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscReal thresh,IS* ISact,IS* ISinact)
63669c03620SShri Abhyankar {
63769c03620SShri Abhyankar   PetscErrorCode ierr;
63869c03620SShri Abhyankar   PetscInt       i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0;
63969c03620SShri Abhyankar   PetscInt       *idx_act,*idx_inact,i1=0,i2=0;
64069c03620SShri Abhyankar   PetscScalar    *db;
64169c03620SShri Abhyankar 
64269c03620SShri Abhyankar   PetscFunctionBegin;
64369c03620SShri Abhyankar 
64469c03620SShri Abhyankar   ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr);
64569c03620SShri Abhyankar   ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr);
64669c03620SShri Abhyankar   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
64769c03620SShri Abhyankar   /* Compute the sizes of the active and inactive sets */
64869c03620SShri Abhyankar   for (i=0; i < nlocal;i++) {
649fca35906SShri Abhyankar     if (PetscAbsScalar(db[i]) <= thresh) nloc_isact++;
650fca35906SShri Abhyankar     else nloc_isinact++;
65169c03620SShri Abhyankar   }
65269c03620SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
65369c03620SShri Abhyankar   ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr);
65469c03620SShri Abhyankar 
65569c03620SShri Abhyankar   /* Creating the indexing arrays */
656fca35906SShri Abhyankar   for(i=0; i < nlocal; i++) {
657fca35906SShri Abhyankar     if (PetscAbsScalar(db[i]) <= thresh) idx_act[i1++] = ilow+i;
658fca35906SShri Abhyankar     else idx_inact[i2++] = ilow+i;
65969c03620SShri Abhyankar   }
66069c03620SShri Abhyankar 
66169c03620SShri Abhyankar   /* Create the index sets */
66269c03620SShri Abhyankar   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr);
66369c03620SShri Abhyankar   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr);
66469c03620SShri Abhyankar 
66569c03620SShri Abhyankar   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
66669c03620SShri Abhyankar   ierr = PetscFree(idx_act);CHKERRQ(ierr);
66769c03620SShri Abhyankar   ierr = PetscFree(idx_inact);CHKERRQ(ierr);
66869c03620SShri Abhyankar   PetscFunctionReturn(0);
66969c03620SShri Abhyankar }
67069c03620SShri Abhyankar 
671dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
672dbd914b8SShri Abhyankar #undef __FUNCT__
673dbd914b8SShri Abhyankar #define __FUNCT__ "SNESVICreateVectors_AS"
674dbd914b8SShri Abhyankar PetscErrorCode SNESVICreateVectors_AS(SNES snes,PetscInt n,Vec* newv)
675dbd914b8SShri Abhyankar {
676dbd914b8SShri Abhyankar   PetscErrorCode ierr;
677dbd914b8SShri Abhyankar   Vec            v;
678dbd914b8SShri Abhyankar 
679dbd914b8SShri Abhyankar   PetscFunctionBegin;
680dbd914b8SShri Abhyankar   ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr);
681dbd914b8SShri Abhyankar   ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr);
682dbd914b8SShri Abhyankar   ierr = VecSetFromOptions(v);CHKERRQ(ierr);
683dbd914b8SShri Abhyankar   *newv = v;
684dbd914b8SShri Abhyankar 
685dbd914b8SShri Abhyankar   PetscFunctionReturn(0);
686dbd914b8SShri Abhyankar }
687dbd914b8SShri Abhyankar 
688dbd914b8SShri Abhyankar 
68969c03620SShri Abhyankar /* Variational Inequality solver using active set method */
69069c03620SShri Abhyankar #undef __FUNCT__
69169c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS"
69269c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes)
69369c03620SShri Abhyankar {
69469c03620SShri Abhyankar   SNES_VI          *vi = (SNES_VI*)snes->data;
69569c03620SShri Abhyankar   PetscErrorCode     ierr;
696408da460SShri Abhyankar   PetscInt           maxits,i,lits,Nis_act=0;
6973b336b49SShri Abhyankar   PetscBool         lssucceed;
69869c03620SShri Abhyankar   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
69969c03620SShri Abhyankar   PetscReal          gnorm,xnorm=0,ynorm;
70069c03620SShri Abhyankar   Vec                Y,X,F,G,W;
70169c03620SShri Abhyankar   KSPConvergedReason kspreason;
70269c03620SShri Abhyankar 
70369c03620SShri Abhyankar   PetscFunctionBegin;
70469c03620SShri Abhyankar   snes->numFailures            = 0;
70569c03620SShri Abhyankar   snes->numLinearSolveFailures = 0;
70669c03620SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
70769c03620SShri Abhyankar 
70869c03620SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
70969c03620SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
71069c03620SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
71169c03620SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
71269c03620SShri Abhyankar   G		= snes->work[1];
71369c03620SShri Abhyankar   W		= snes->work[2];
71469c03620SShri Abhyankar 
71569c03620SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
71669c03620SShri Abhyankar   snes->iter = 0;
71769c03620SShri Abhyankar   snes->norm = 0.0;
71869c03620SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
71969c03620SShri Abhyankar 
72069c03620SShri Abhyankar   ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr);
72169c03620SShri Abhyankar   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
72269c03620SShri Abhyankar   if (snes->domainerror) {
72369c03620SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
72469c03620SShri Abhyankar     PetscFunctionReturn(0);
72569c03620SShri Abhyankar   }
72669c03620SShri Abhyankar    /* Compute Merit function */
72769c03620SShri Abhyankar   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
72869c03620SShri Abhyankar 
72969c03620SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
73069c03620SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
73169c03620SShri Abhyankar   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
73269c03620SShri Abhyankar 
73369c03620SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
73469c03620SShri Abhyankar   snes->norm = vi->phinorm;
73569c03620SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
73669c03620SShri Abhyankar   SNESLogConvHistory(snes,vi->phinorm,0);
73769c03620SShri Abhyankar   SNESMonitor(snes,0,vi->phinorm);
73869c03620SShri Abhyankar 
73969c03620SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
74069c03620SShri Abhyankar   snes->ttol = vi->phinorm*snes->rtol;
74169c03620SShri Abhyankar   /* test convergence */
74269c03620SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
74369c03620SShri Abhyankar   if (snes->reason) PetscFunctionReturn(0);
74469c03620SShri Abhyankar 
74569c03620SShri Abhyankar   for (i=0; i<maxits; i++) {
74669c03620SShri Abhyankar 
747a79edbefSShri Abhyankar     IS                 IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
748a79edbefSShri Abhyankar     PetscReal          thresh,J_norm1;
749dbd914b8SShri Abhyankar     VecScatter         scat_act,scat_inact;
750408da460SShri Abhyankar     PetscInt           nis_act,nis_inact,Nis_act_prev;
751dbd914b8SShri Abhyankar     Vec                Da_act,Da_inact,Db_inact;
752dbd914b8SShri Abhyankar     Vec                Y_act,Y_inact,phi_act,phi_inact;
753d76c674bSShri Abhyankar     Mat                jac_inact_inact,jac_inact_act,prejac_inact_inact;
754a79edbefSShri Abhyankar 
75569c03620SShri Abhyankar     /* Call general purpose update function */
75669c03620SShri Abhyankar     if (snes->ops->update) {
75769c03620SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
75869c03620SShri Abhyankar     }
75969c03620SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
76069c03620SShri Abhyankar     /* Compute the threshold value for creating active and inactive sets */
76169c03620SShri Abhyankar     ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr);
76269c03620SShri Abhyankar     thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1);
763a79edbefSShri Abhyankar 
764a79edbefSShri Abhyankar     /* Compute B-subdifferential vectors Da and Db */
765a79edbefSShri Abhyankar     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
766dbd914b8SShri Abhyankar 
76769c03620SShri Abhyankar     /* Create active and inactive index sets */
76869c03620SShri Abhyankar     ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr);
76969c03620SShri Abhyankar 
770408da460SShri Abhyankar     Nis_act_prev = Nis_act;
771408da460SShri Abhyankar     /* Get sizes of active and inactive sets */
772dbd914b8SShri Abhyankar     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
773dbd914b8SShri Abhyankar     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
774408da460SShri Abhyankar     ierr = ISGetSize(IS_act,&Nis_act);CHKERRQ(ierr);
775dbd914b8SShri Abhyankar 
776d76c674bSShri Abhyankar     ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr);
777d76c674bSShri Abhyankar 
778dbd914b8SShri Abhyankar     /* Create active and inactive set vectors */
779dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_act,&Da_act);CHKERRQ(ierr);
780dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_inact,&Da_inact);CHKERRQ(ierr);
781dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_inact,&Db_inact);CHKERRQ(ierr);
782dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_act,&phi_act);CHKERRQ(ierr);
783dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr);
784dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr);
785dbd914b8SShri Abhyankar     ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
786dbd914b8SShri Abhyankar 
787dbd914b8SShri Abhyankar     /* Create inactive set submatrices */
788dbd914b8SShri Abhyankar     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_act,MAT_INITIAL_MATRIX,&jac_inact_act);CHKERRQ(ierr);
789dbd914b8SShri Abhyankar     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
790dbd914b8SShri Abhyankar 
791dbd914b8SShri Abhyankar     /* Create scatter contexts */
792dbd914b8SShri Abhyankar     ierr = VecScatterCreate(vi->Da,IS_act,Da_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
793dbd914b8SShri Abhyankar     ierr = VecScatterCreate(vi->Da,IS_inact,Da_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
794dbd914b8SShri Abhyankar 
795dbd914b8SShri Abhyankar     /* Do a vec scatter to active and inactive set vectors */
796dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
797dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
798dbd914b8SShri Abhyankar 
799dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
800dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
801dbd914b8SShri Abhyankar 
802dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
803dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
804dbd914b8SShri Abhyankar 
805dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
806dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
807dbd914b8SShri Abhyankar 
808dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
809dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
810dbd914b8SShri Abhyankar 
811dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
812dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
813dbd914b8SShri Abhyankar 
814dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
815dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
816dbd914b8SShri Abhyankar 
817dbd914b8SShri Abhyankar     /* Active set direction */
818dbd914b8SShri Abhyankar     ierr = VecPointwiseDivide(Y_act,phi_act,Da_act);CHKERRQ(ierr);
819d76c674bSShri Abhyankar     /* inactive set jacobian and preconditioner */
820dbd914b8SShri Abhyankar     ierr = VecPointwiseDivide(Da_inact,Da_inact,Db_inact);CHKERRQ(ierr);
821dbd914b8SShri Abhyankar     ierr = MatDiagonalSet(jac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr);
822d76c674bSShri Abhyankar     if (snes->jacobian != snes->jacobian_pre) {
823d76c674bSShri Abhyankar       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
824d76c674bSShri Abhyankar       ierr = MatDiagonalSet(prejac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr);
825d76c674bSShri Abhyankar     } else prejac_inact_inact = jac_inact_inact;
826d76c674bSShri Abhyankar 
827dbd914b8SShri Abhyankar     /* right hand side */
828dbd914b8SShri Abhyankar     ierr = VecPointwiseDivide(phi_inact,phi_inact,Db_inact);CHKERRQ(ierr);
829dbd914b8SShri Abhyankar     ierr = MatMult(jac_inact_act,Y_act,Db_inact);CHKERRQ(ierr);
830dbd914b8SShri Abhyankar     ierr = VecAXPY(phi_inact,-1.0,Db_inact);CHKERRQ(ierr);
831dbd914b8SShri Abhyankar 
832408da460SShri Abhyankar     if ((i != 0) && (Nis_act != Nis_act_prev)) {
833408da460SShri Abhyankar       KSP kspnew,snesksp;
834408da460SShri Abhyankar       PC  pcnew;
835408da460SShri Abhyankar       /* The active and inactive set sizes have changed so need to create a new snes->ksp object */
836408da460SShri Abhyankar       ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
837408da460SShri Abhyankar       ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr);
838408da460SShri Abhyankar       /* Copy over snes->ksp info */
839408da460SShri Abhyankar       kspnew->pc_side = snesksp->pc_side;
840408da460SShri Abhyankar       kspnew->rtol    = snesksp->rtol;
841408da460SShri Abhyankar       kspnew->abstol    = snesksp->abstol;
842408da460SShri Abhyankar       kspnew->max_it  = snesksp->max_it;
843408da460SShri Abhyankar       ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr);
844408da460SShri Abhyankar       ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr);
845408da460SShri Abhyankar       ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr);
846408da460SShri Abhyankar       ierr = KSPDestroy(snesksp);CHKERRQ(ierr);
847408da460SShri Abhyankar       snes->ksp = kspnew;
848408da460SShri Abhyankar       ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr);
849408da460SShri Abhyankar       ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);
850408da460SShri Abhyankar     }
851408da460SShri Abhyankar 
852d76c674bSShri Abhyankar     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
853dbd914b8SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr);
854dbd914b8SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
8553b336b49SShri Abhyankar     if (kspreason < 0) {
8563b336b49SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
8573b336b49SShri 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);
8583b336b49SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
8593b336b49SShri Abhyankar         break;
8603b336b49SShri Abhyankar       }
8613b336b49SShri Abhyankar      }
862dbd914b8SShri Abhyankar 
863dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
864dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
865dbd914b8SShri Abhyankar     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
866dbd914b8SShri Abhyankar     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
867dbd914b8SShri Abhyankar 
868dbd914b8SShri Abhyankar     ierr = VecDestroy(Da_act);CHKERRQ(ierr);
869dbd914b8SShri Abhyankar     ierr = VecDestroy(Da_inact);CHKERRQ(ierr);
870dbd914b8SShri Abhyankar     ierr = VecDestroy(Db_inact);CHKERRQ(ierr);
871dbd914b8SShri Abhyankar     ierr = VecDestroy(phi_act);CHKERRQ(ierr);
872dbd914b8SShri Abhyankar     ierr = VecDestroy(phi_inact);CHKERRQ(ierr);
873dbd914b8SShri Abhyankar     ierr = VecDestroy(Y_act);CHKERRQ(ierr);
874dbd914b8SShri Abhyankar     ierr = VecDestroy(Y_inact);CHKERRQ(ierr);
875dbd914b8SShri Abhyankar     ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr);
876dbd914b8SShri Abhyankar     ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr);
877730c24dcSShri Abhyankar     ierr = ISDestroy(IS_act);CHKERRQ(ierr);
878730c24dcSShri Abhyankar     ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
879bac2f86eSShri Abhyankar     ierr = MatDestroy(jac_inact_act);CHKERRQ(ierr);
880bac2f86eSShri Abhyankar     ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
881d76c674bSShri Abhyankar     if (snes->jacobian != snes->jacobian_pre) {
882d76c674bSShri Abhyankar       ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr);
883d76c674bSShri Abhyankar     }
884730c24dcSShri Abhyankar 
88569c03620SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
88669c03620SShri Abhyankar     snes->linear_its += lits;
88769c03620SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
88869c03620SShri Abhyankar     /*
88969c03620SShri Abhyankar     if (vi->precheckstep) {
89069c03620SShri Abhyankar       PetscBool changed_y = PETSC_FALSE;
89169c03620SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
89269c03620SShri Abhyankar     }
89369c03620SShri Abhyankar 
89469c03620SShri Abhyankar     if (PetscLogPrintInfo){
89569c03620SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
89669c03620SShri Abhyankar     }
89769c03620SShri Abhyankar     */
89869c03620SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
89969c03620SShri Abhyankar          Y <- X - lambda*Y
90069c03620SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
90169c03620SShri Abhyankar     */
90269c03620SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
90369c03620SShri Abhyankar     ynorm = 1; gnorm = vi->phinorm;
9043b336b49SShri Abhyankar     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
90569c03620SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
90669c03620SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
90769c03620SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
90869c03620SShri Abhyankar     if (snes->domainerror) {
90969c03620SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
91069c03620SShri Abhyankar       PetscFunctionReturn(0);
91169c03620SShri Abhyankar     }
91269c03620SShri Abhyankar     if (!lssucceed) {
91369c03620SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
91469c03620SShri Abhyankar 	PetscBool ismin;
91569c03620SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
91669c03620SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
91769c03620SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
91869c03620SShri Abhyankar         break;
91969c03620SShri Abhyankar       }
92069c03620SShri Abhyankar     }
92169c03620SShri Abhyankar     /* Update function and solution vectors */
92269c03620SShri Abhyankar     vi->phinorm = gnorm;
92369c03620SShri Abhyankar     vi->merit = 0.5*vi->phinorm*vi->phinorm;
92469c03620SShri Abhyankar     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
92569c03620SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
92669c03620SShri Abhyankar     /* Monitor convergence */
92769c03620SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
92869c03620SShri Abhyankar     snes->iter = i+1;
92969c03620SShri Abhyankar     snes->norm = vi->phinorm;
93069c03620SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
93169c03620SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
93269c03620SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
93369c03620SShri Abhyankar     /* Test for convergence, xnorm = || X || */
93469c03620SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
93569c03620SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
93669c03620SShri Abhyankar     if (snes->reason) break;
93769c03620SShri Abhyankar   }
93869c03620SShri Abhyankar   if (i == maxits) {
93969c03620SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
94069c03620SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
94169c03620SShri Abhyankar   }
94269c03620SShri Abhyankar   PetscFunctionReturn(0);
94369c03620SShri Abhyankar }
94469c03620SShri Abhyankar 
9452b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
9462b4ed282SShri Abhyankar /*
9472b4ed282SShri Abhyankar    SNESSetUp_VI - Sets up the internal data structures for the later use
9482b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
9492b4ed282SShri Abhyankar 
9502b4ed282SShri Abhyankar    Input Parameter:
9512b4ed282SShri Abhyankar .  snes - the SNES context
9522b4ed282SShri Abhyankar .  x - the solution vector
9532b4ed282SShri Abhyankar 
9542b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
9552b4ed282SShri Abhyankar 
9562b4ed282SShri Abhyankar    Notes:
9572b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
9582b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
9592b4ed282SShri Abhyankar    the call to SNESSolve().
9602b4ed282SShri Abhyankar  */
9612b4ed282SShri Abhyankar #undef __FUNCT__
9622b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
9632b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
9642b4ed282SShri Abhyankar {
9652b4ed282SShri Abhyankar   PetscErrorCode ierr;
9662b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*) snes->data;
9672b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
9682b4ed282SShri Abhyankar 
9692b4ed282SShri Abhyankar   PetscFunctionBegin;
9702b4ed282SShri Abhyankar   if (!snes->vec_sol_update) {
9712b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
9722b4ed282SShri Abhyankar     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
9732b4ed282SShri Abhyankar   }
9742b4ed282SShri Abhyankar   if (!snes->work) {
9752b4ed282SShri Abhyankar     snes->nwork = 3;
9762b4ed282SShri Abhyankar     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
9772b4ed282SShri Abhyankar     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
9782b4ed282SShri Abhyankar   }
9792b4ed282SShri Abhyankar 
9802b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr);
9812b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr);
9822b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr);
9832b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
9842b4ed282SShri Abhyankar   ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr);
9852b4ed282SShri Abhyankar 
9862b4ed282SShri Abhyankar   /* If the lower and upper bound on variables are not set, set it to
9872b4ed282SShri Abhyankar      -Inf and Inf */
9882b4ed282SShri Abhyankar   if (!vi->xl && !vi->xu) {
9892b4ed282SShri Abhyankar     vi->usersetxbounds = PETSC_FALSE;
9902b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr);
9912b4ed282SShri Abhyankar     ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr);
9922b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr);
9932b4ed282SShri Abhyankar     ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr);
9942b4ed282SShri Abhyankar   } else {
9952b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
9962b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
9972b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
9982b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
9992b4ed282SShri 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]))
10002b4ed282SShri 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.");
10012b4ed282SShri Abhyankar   }
10022b4ed282SShri Abhyankar 
10032b4ed282SShri Abhyankar   vi->computeuserfunction = snes->ops->computefunction;
10041f79c099SShri Abhyankar   snes->ops->computefunction = SNESVIComputeFunction;
1005a79edbefSShri Abhyankar 
10062b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10072b4ed282SShri Abhyankar }
10082b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10092b4ed282SShri Abhyankar /*
10102b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
10112b4ed282SShri Abhyankar    with SNESCreate_VI().
10122b4ed282SShri Abhyankar 
10132b4ed282SShri Abhyankar    Input Parameter:
10142b4ed282SShri Abhyankar .  snes - the SNES context
10152b4ed282SShri Abhyankar 
10162b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
10172b4ed282SShri Abhyankar  */
10182b4ed282SShri Abhyankar #undef __FUNCT__
10192b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
10202b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
10212b4ed282SShri Abhyankar {
10222b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*) snes->data;
10232b4ed282SShri Abhyankar   PetscErrorCode ierr;
10242b4ed282SShri Abhyankar 
10252b4ed282SShri Abhyankar   PetscFunctionBegin;
10262b4ed282SShri Abhyankar   if (snes->vec_sol_update) {
10272b4ed282SShri Abhyankar     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
10282b4ed282SShri Abhyankar     snes->vec_sol_update = PETSC_NULL;
10292b4ed282SShri Abhyankar   }
10302b4ed282SShri Abhyankar   if (snes->nwork) {
10312b4ed282SShri Abhyankar     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
10322b4ed282SShri Abhyankar     snes->nwork = 0;
10332b4ed282SShri Abhyankar     snes->work  = PETSC_NULL;
10342b4ed282SShri Abhyankar   }
10352b4ed282SShri Abhyankar 
10362b4ed282SShri Abhyankar   /* clear vectors */
10372b4ed282SShri Abhyankar   ierr = VecDestroy(vi->phi); CHKERRQ(ierr);
10382b4ed282SShri Abhyankar   ierr = VecDestroy(vi->Da); CHKERRQ(ierr);
10392b4ed282SShri Abhyankar   ierr = VecDestroy(vi->Db); CHKERRQ(ierr);
10402b4ed282SShri Abhyankar   ierr = VecDestroy(vi->z); CHKERRQ(ierr);
10412b4ed282SShri Abhyankar   ierr = VecDestroy(vi->t); CHKERRQ(ierr);
10422b4ed282SShri Abhyankar   if (!vi->usersetxbounds) {
10432b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xl); CHKERRQ(ierr);
10442b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xu); CHKERRQ(ierr);
10452b4ed282SShri Abhyankar   }
1046c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1047c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
10482b4ed282SShri Abhyankar   }
10492b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
10502b4ed282SShri Abhyankar 
10512b4ed282SShri Abhyankar   /* clear composed functions */
10522b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
1053c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
10542b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10552b4ed282SShri Abhyankar }
10562b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10572b4ed282SShri Abhyankar #undef __FUNCT__
10582b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
10592b4ed282SShri Abhyankar 
10602b4ed282SShri Abhyankar /*
10612b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c
10622b4ed282SShri Abhyankar 
10632b4ed282SShri Abhyankar */
1064ace3abfcSBarry 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)
10652b4ed282SShri Abhyankar {
10662b4ed282SShri Abhyankar   PetscErrorCode ierr;
10672b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1068ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
10692b4ed282SShri Abhyankar 
10702b4ed282SShri Abhyankar   PetscFunctionBegin;
10712b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
10722b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10732b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
10742b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
10752b4ed282SShri Abhyankar   if (vi->postcheckstep) {
10762b4ed282SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
10772b4ed282SShri Abhyankar   }
10782b4ed282SShri Abhyankar   if (changed_y) {
10792b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
10802b4ed282SShri Abhyankar   }
10812b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
10822b4ed282SShri Abhyankar   if (!snes->domainerror) {
10832b4ed282SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
10842b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
10852b4ed282SShri Abhyankar   }
1086c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1087c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1088c92abb8aSShri Abhyankar   }
10892b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10902b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10912b4ed282SShri Abhyankar }
10922b4ed282SShri Abhyankar 
10932b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10942b4ed282SShri Abhyankar #undef __FUNCT__
10952b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
10962b4ed282SShri Abhyankar 
10972b4ed282SShri Abhyankar /*
10982b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
10992b4ed282SShri Abhyankar */
1100ace3abfcSBarry 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)
11012b4ed282SShri Abhyankar {
11022b4ed282SShri Abhyankar   PetscErrorCode ierr;
11032b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1104ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
11052b4ed282SShri Abhyankar 
11062b4ed282SShri Abhyankar   PetscFunctionBegin;
11072b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
11082b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11092b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
11102b4ed282SShri Abhyankar   if (vi->postcheckstep) {
11112b4ed282SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
11122b4ed282SShri Abhyankar   }
11132b4ed282SShri Abhyankar   if (changed_y) {
11142b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
11152b4ed282SShri Abhyankar   }
11162b4ed282SShri Abhyankar 
11172b4ed282SShri Abhyankar   /* don't evaluate function the last time through */
11182b4ed282SShri Abhyankar   if (snes->iter < snes->max_its-1) {
11192b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
11202b4ed282SShri Abhyankar   }
11212b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11222b4ed282SShri Abhyankar   PetscFunctionReturn(0);
11232b4ed282SShri Abhyankar }
11242b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
11252b4ed282SShri Abhyankar #undef __FUNCT__
11262b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
11272b4ed282SShri Abhyankar /*
11282b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c
11292b4ed282SShri Abhyankar */
1130ace3abfcSBarry 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)
11312b4ed282SShri Abhyankar {
11322b4ed282SShri Abhyankar   /*
11332b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
11342b4ed282SShri Abhyankar      minimization problem:
11352b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
11362b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
11372b4ed282SShri Abhyankar      This function z(x) is same as the merit function for the complementarity problem.
11382b4ed282SShri Abhyankar    */
11392b4ed282SShri Abhyankar 
11402b4ed282SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
11412b4ed282SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
11422b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11432b4ed282SShri Abhyankar   PetscScalar    cinitslope;
11442b4ed282SShri Abhyankar #endif
11452b4ed282SShri Abhyankar   PetscErrorCode ierr;
11462b4ed282SShri Abhyankar   PetscInt       count;
11472b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*)snes->data;
1148ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
11492b4ed282SShri Abhyankar   MPI_Comm       comm;
11502b4ed282SShri Abhyankar 
11512b4ed282SShri Abhyankar   PetscFunctionBegin;
11522b4ed282SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
11532b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11542b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
11552b4ed282SShri Abhyankar 
11562b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
11572b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
1158c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1159c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
11602b4ed282SShri Abhyankar     }
11612b4ed282SShri Abhyankar     *gnorm = fnorm;
11622b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
11632b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
11642b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
11652b4ed282SShri Abhyankar     goto theend1;
11662b4ed282SShri Abhyankar   }
11672b4ed282SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1168c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1169c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
11702b4ed282SShri Abhyankar     }
11712b4ed282SShri Abhyankar     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
11722b4ed282SShri Abhyankar     *ynorm = vi->maxstep;
11732b4ed282SShri Abhyankar   }
11742b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
11752b4ed282SShri Abhyankar   minlambda = vi->minlambda/rellength;
11763b336b49SShri Abhyankar   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
11772b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11783b336b49SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
11792b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
11802b4ed282SShri Abhyankar #else
11813b336b49SShri Abhyankar   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
11822b4ed282SShri Abhyankar #endif
11832b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
11842b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
11852b4ed282SShri Abhyankar 
11862b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
11872b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
11882b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
11892b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
11902b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
11912b4ed282SShri Abhyankar     goto theend1;
11922b4ed282SShri Abhyankar   }
11932b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
11942b4ed282SShri Abhyankar   if (snes->domainerror) {
11952b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11962b4ed282SShri Abhyankar     PetscFunctionReturn(0);
11972b4ed282SShri Abhyankar   }
11982b4ed282SShri Abhyankar   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
11992b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
12002b4ed282SShri Abhyankar   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
12012b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1202c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1203c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
12042b4ed282SShri Abhyankar     }
12052b4ed282SShri Abhyankar     goto theend1;
12062b4ed282SShri Abhyankar   }
12072b4ed282SShri Abhyankar 
12082b4ed282SShri Abhyankar   /* Fit points with quadratic */
12092b4ed282SShri Abhyankar   lambda     = 1.0;
12102b4ed282SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
12112b4ed282SShri Abhyankar   lambdaprev = lambda;
12122b4ed282SShri Abhyankar   gnormprev  = *gnorm;
12132b4ed282SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
12142b4ed282SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
12152b4ed282SShri Abhyankar   else                         lambda = lambdatemp;
12162b4ed282SShri Abhyankar 
12172b4ed282SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
12182b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
12192b4ed282SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
12202b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
12212b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
12222b4ed282SShri Abhyankar     goto theend1;
12232b4ed282SShri Abhyankar   }
12242b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12252b4ed282SShri Abhyankar   if (snes->domainerror) {
12262b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
12272b4ed282SShri Abhyankar     PetscFunctionReturn(0);
12282b4ed282SShri Abhyankar   }
12292b4ed282SShri Abhyankar   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12302b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1231c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1232c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
12332b4ed282SShri Abhyankar   }
12342b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1235c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1236c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
12372b4ed282SShri Abhyankar     }
12382b4ed282SShri Abhyankar     goto theend1;
12392b4ed282SShri Abhyankar   }
12402b4ed282SShri Abhyankar 
12412b4ed282SShri Abhyankar   /* Fit points with cubic */
12422b4ed282SShri Abhyankar   count = 1;
12432b4ed282SShri Abhyankar   while (PETSC_TRUE) {
12442b4ed282SShri Abhyankar     if (lambda <= minlambda) {
1245c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1246c92abb8aSShri Abhyankar 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
1247c92abb8aSShri 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);
12482b4ed282SShri Abhyankar       }
12492b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
12502b4ed282SShri Abhyankar       break;
12512b4ed282SShri Abhyankar     }
12522b4ed282SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
12532b4ed282SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
12542b4ed282SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
12552b4ed282SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
12562b4ed282SShri Abhyankar     d  = b*b - 3*a*initslope;
12572b4ed282SShri Abhyankar     if (d < 0.0) d = 0.0;
12582b4ed282SShri Abhyankar     if (a == 0.0) {
12592b4ed282SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
12602b4ed282SShri Abhyankar     } else {
12612b4ed282SShri Abhyankar       lambdatemp = (-b + sqrt(d))/(3.0*a);
12622b4ed282SShri Abhyankar     }
12632b4ed282SShri Abhyankar     lambdaprev = lambda;
12642b4ed282SShri Abhyankar     gnormprev  = *gnorm;
12652b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
12662b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
12672b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
12682b4ed282SShri Abhyankar     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
12692b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
12702b4ed282SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
12712b4ed282SShri 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);
12722b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
12732b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
12742b4ed282SShri Abhyankar       break;
12752b4ed282SShri Abhyankar     }
12762b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12772b4ed282SShri Abhyankar     if (snes->domainerror) {
12782b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
12792b4ed282SShri Abhyankar       PetscFunctionReturn(0);
12802b4ed282SShri Abhyankar     }
12812b4ed282SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12822b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
12832b4ed282SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
1284c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
12852b4ed282SShri Abhyankar 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
12862b4ed282SShri Abhyankar       }
12872b4ed282SShri Abhyankar       break;
12882b4ed282SShri Abhyankar     } else {
1289c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
12902b4ed282SShri Abhyankar         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
12912b4ed282SShri Abhyankar       }
12922b4ed282SShri Abhyankar     }
12932b4ed282SShri Abhyankar     count++;
12942b4ed282SShri Abhyankar   }
12952b4ed282SShri Abhyankar   theend1:
12962b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
12972b4ed282SShri Abhyankar   if (vi->postcheckstep && *flag) {
12982b4ed282SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
12992b4ed282SShri Abhyankar     if (changed_y) {
13002b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
13012b4ed282SShri Abhyankar     }
13022b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
13032b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13042b4ed282SShri Abhyankar       if (snes->domainerror) {
13052b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13062b4ed282SShri Abhyankar         PetscFunctionReturn(0);
13072b4ed282SShri Abhyankar       }
13082b4ed282SShri Abhyankar       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
13092b4ed282SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
13102b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
13112b4ed282SShri Abhyankar       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
13122b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
13132b4ed282SShri Abhyankar     }
13142b4ed282SShri Abhyankar   }
13152b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13162b4ed282SShri Abhyankar   PetscFunctionReturn(0);
13172b4ed282SShri Abhyankar }
13182b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
13192b4ed282SShri Abhyankar #undef __FUNCT__
13202b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
13212b4ed282SShri Abhyankar /*
13222b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c
13232b4ed282SShri Abhyankar */
1324ace3abfcSBarry 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)
13252b4ed282SShri Abhyankar {
13262b4ed282SShri Abhyankar   /*
13272b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
13282b4ed282SShri Abhyankar      minimization problem:
13292b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
13302b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
13312b4ed282SShri Abhyankar      z(x) is the same as the merit function for the complementarity problem
13322b4ed282SShri Abhyankar    */
13332b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
13342b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
13352b4ed282SShri Abhyankar   PetscScalar    cinitslope;
13362b4ed282SShri Abhyankar #endif
13372b4ed282SShri Abhyankar   PetscErrorCode ierr;
13382b4ed282SShri Abhyankar   PetscInt       count;
13392b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1340ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
13412b4ed282SShri Abhyankar 
13422b4ed282SShri Abhyankar   PetscFunctionBegin;
13432b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13442b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
13452b4ed282SShri Abhyankar 
13462b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
13472b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
1348c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1349c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
13502b4ed282SShri Abhyankar     }
13512b4ed282SShri Abhyankar     *gnorm = fnorm;
13522b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
13532b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
13542b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
13552b4ed282SShri Abhyankar     goto theend2;
13562b4ed282SShri Abhyankar   }
13572b4ed282SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
13582b4ed282SShri Abhyankar     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
13592b4ed282SShri Abhyankar     *ynorm = vi->maxstep;
13602b4ed282SShri Abhyankar   }
13612b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
13622b4ed282SShri Abhyankar   minlambda = vi->minlambda/rellength;
13633b336b49SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
13642b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
13652b4ed282SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
13662b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
13672b4ed282SShri Abhyankar #else
13683b336b49SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
13692b4ed282SShri Abhyankar #endif
13702b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
13712b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
13722b4ed282SShri Abhyankar   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
13732b4ed282SShri Abhyankar 
13742b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
13752b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
13762b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
13772b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
13782b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
13792b4ed282SShri Abhyankar     goto theend2;
13802b4ed282SShri Abhyankar   }
13812b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
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 = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13872b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
13882b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1389c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1390c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
13912b4ed282SShri Abhyankar     }
13922b4ed282SShri Abhyankar     goto theend2;
13932b4ed282SShri Abhyankar   }
13942b4ed282SShri Abhyankar 
13952b4ed282SShri Abhyankar   /* Fit points with quadratic */
13962b4ed282SShri Abhyankar   lambda = 1.0;
13972b4ed282SShri Abhyankar   count = 1;
13982b4ed282SShri Abhyankar   while (PETSC_TRUE) {
13992b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
1400c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1401c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
1402c92abb8aSShri 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);
14032b4ed282SShri Abhyankar       }
14042b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
14052b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
14062b4ed282SShri Abhyankar       break;
14072b4ed282SShri Abhyankar     }
14082b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
14092b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
14102b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
14112b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
14122b4ed282SShri Abhyankar 
14132b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
14142b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
14152b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
14162b4ed282SShri 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);
14172b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
14182b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
14192b4ed282SShri Abhyankar       break;
14202b4ed282SShri Abhyankar     }
14212b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
14222b4ed282SShri Abhyankar     if (snes->domainerror) {
14232b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14242b4ed282SShri Abhyankar       PetscFunctionReturn(0);
14252b4ed282SShri Abhyankar     }
14262b4ed282SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
14272b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
14282b4ed282SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1429c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1430c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
14312b4ed282SShri Abhyankar       }
14322b4ed282SShri Abhyankar       break;
14332b4ed282SShri Abhyankar     }
14342b4ed282SShri Abhyankar     count++;
14352b4ed282SShri Abhyankar   }
14362b4ed282SShri Abhyankar   theend2:
14372b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
14382b4ed282SShri Abhyankar   if (vi->postcheckstep) {
14392b4ed282SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
14402b4ed282SShri Abhyankar     if (changed_y) {
14412b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
14422b4ed282SShri Abhyankar     }
14432b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
14442b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
14452b4ed282SShri Abhyankar       if (snes->domainerror) {
14462b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14472b4ed282SShri Abhyankar         PetscFunctionReturn(0);
14482b4ed282SShri Abhyankar       }
14492b4ed282SShri Abhyankar       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
14502b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
14512b4ed282SShri Abhyankar       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
14522b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
14532b4ed282SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
14542b4ed282SShri Abhyankar     }
14552b4ed282SShri Abhyankar   }
14562b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14572b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14582b4ed282SShri Abhyankar }
14592b4ed282SShri Abhyankar 
1460ace3abfcSBarry 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*/
14612b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
14622b4ed282SShri Abhyankar EXTERN_C_BEGIN
14632b4ed282SShri Abhyankar #undef __FUNCT__
14642b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI"
14652b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
14662b4ed282SShri Abhyankar {
14672b4ed282SShri Abhyankar   PetscFunctionBegin;
14682b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->LineSearch = func;
14692b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->lsP        = lsctx;
14702b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14712b4ed282SShri Abhyankar }
14722b4ed282SShri Abhyankar EXTERN_C_END
14732b4ed282SShri Abhyankar 
14742b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
14752b4ed282SShri Abhyankar EXTERN_C_BEGIN
14762b4ed282SShri Abhyankar #undef __FUNCT__
14772b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
1478ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
14792b4ed282SShri Abhyankar {
14802b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
14812b4ed282SShri Abhyankar   PetscErrorCode ierr;
14822b4ed282SShri Abhyankar 
14832b4ed282SShri Abhyankar   PetscFunctionBegin;
1484c92abb8aSShri Abhyankar   if (flg && !vi->lsmonitor) {
1485c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr);
1486c92abb8aSShri Abhyankar   } else if (!flg && vi->lsmonitor) {
1487c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
14882b4ed282SShri Abhyankar   }
14892b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14902b4ed282SShri Abhyankar }
14912b4ed282SShri Abhyankar EXTERN_C_END
14922b4ed282SShri Abhyankar 
14932b4ed282SShri Abhyankar /*
14942b4ed282SShri Abhyankar    SNESView_VI - Prints info from the SNESVI data structure.
14952b4ed282SShri Abhyankar 
14962b4ed282SShri Abhyankar    Input Parameters:
14972b4ed282SShri Abhyankar .  SNES - the SNES context
14982b4ed282SShri Abhyankar .  viewer - visualization context
14992b4ed282SShri Abhyankar 
15002b4ed282SShri Abhyankar    Application Interface Routine: SNESView()
15012b4ed282SShri Abhyankar */
15022b4ed282SShri Abhyankar #undef __FUNCT__
15032b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI"
15042b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
15052b4ed282SShri Abhyankar {
15062b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
150778c4b9d3SShri Abhyankar   const char     *cstr,*tstr;
15082b4ed282SShri Abhyankar   PetscErrorCode ierr;
1509ace3abfcSBarry Smith   PetscBool     iascii;
15102b4ed282SShri Abhyankar 
15112b4ed282SShri Abhyankar   PetscFunctionBegin;
15122b4ed282SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
15132b4ed282SShri Abhyankar   if (iascii) {
15142b4ed282SShri Abhyankar     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
15152b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
15162b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
15172b4ed282SShri Abhyankar     else                                                cstr = "unknown";
151878c4b9d3SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_SS)      tstr = "Semismooth";
151978c4b9d3SShri Abhyankar     else if (snes->ops->solve == SNESSolveVI_AS)  tstr = "Active Set";
152078c4b9d3SShri Abhyankar     else                                         tstr = "unknown";
152178c4b9d3SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  VI algorithm: %s\n",tstr);CHKERRQ(ierr);
15222b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
15232b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
15242b4ed282SShri Abhyankar   } else {
15252b4ed282SShri Abhyankar     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
15262b4ed282SShri Abhyankar   }
15272b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15282b4ed282SShri Abhyankar }
15292b4ed282SShri Abhyankar 
15302b4ed282SShri Abhyankar /*
15312b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
15322b4ed282SShri Abhyankar 
15332b4ed282SShri Abhyankar    Input Parameters:
15342b4ed282SShri Abhyankar .  snes - the SNES context.
15352b4ed282SShri Abhyankar .  xl   - lower bound.
15362b4ed282SShri Abhyankar .  xu   - upper bound.
15372b4ed282SShri Abhyankar 
15382b4ed282SShri Abhyankar    Notes:
15392b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
1540*84c105d7SBarry Smith    PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp().
1541*84c105d7SBarry Smith 
15422b4ed282SShri Abhyankar */
15432b4ed282SShri Abhyankar 
15442b4ed282SShri Abhyankar #undef __FUNCT__
15452b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds"
154669c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
15472b4ed282SShri Abhyankar {
15482b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
15492b4ed282SShri Abhyankar 
15502b4ed282SShri Abhyankar   PetscFunctionBegin;
15512b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
15522b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
15532b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
15542b4ed282SShri Abhyankar 
15552b4ed282SShri Abhyankar   /* Check if SNESSetFunction is called */
15562b4ed282SShri Abhyankar   if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
15572b4ed282SShri Abhyankar 
15582b4ed282SShri Abhyankar   /* Check if the vector sizes are compatible for lower and upper bounds */
15592b4ed282SShri 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);
15602b4ed282SShri 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);
15612b4ed282SShri Abhyankar   vi->usersetxbounds = PETSC_TRUE;
15622b4ed282SShri Abhyankar   vi->xl = xl;
15632b4ed282SShri Abhyankar   vi->xu = xu;
15642b4ed282SShri Abhyankar 
15652b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15662b4ed282SShri Abhyankar }
15672b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15682b4ed282SShri Abhyankar /*
15692b4ed282SShri Abhyankar    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
15702b4ed282SShri Abhyankar 
15712b4ed282SShri Abhyankar    Input Parameter:
15722b4ed282SShri Abhyankar .  snes - the SNES context
15732b4ed282SShri Abhyankar 
15742b4ed282SShri Abhyankar    Application Interface Routine: SNESSetFromOptions()
15752b4ed282SShri Abhyankar */
15762b4ed282SShri Abhyankar #undef __FUNCT__
15772b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI"
15782b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
15792b4ed282SShri Abhyankar {
15802b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
15812b4ed282SShri Abhyankar   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
158269c03620SShri Abhyankar   const char     *vies[] = {"ss","as"};
15832b4ed282SShri Abhyankar   PetscErrorCode ierr;
15842b4ed282SShri Abhyankar   PetscInt       indx;
158569c03620SShri Abhyankar   PetscBool     flg,set,flg2;
15862b4ed282SShri Abhyankar 
15872b4ed282SShri Abhyankar   PetscFunctionBegin;
15882b4ed282SShri Abhyankar     ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
15892b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
15902b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
15912b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
15922b4ed282SShri Abhyankar     ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
1593acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
15942b4ed282SShri Abhyankar     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
159569c03620SShri Abhyankar     ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr);
159669c03620SShri Abhyankar     if (flg2) {
159769c03620SShri Abhyankar       switch (indx) {
159869c03620SShri Abhyankar       case 0:
159969c03620SShri Abhyankar 	snes->ops->solve = SNESSolveVI_SS;
160069c03620SShri Abhyankar 	break;
160169c03620SShri Abhyankar       case 1:
160269c03620SShri Abhyankar 	snes->ops->solve = SNESSolveVI_AS;
160369c03620SShri Abhyankar 	break;
160469c03620SShri Abhyankar       }
160569c03620SShri Abhyankar     }
1606c92abb8aSShri Abhyankar     ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
16072b4ed282SShri Abhyankar     if (flg) {
16082b4ed282SShri Abhyankar       switch (indx) {
16092b4ed282SShri Abhyankar       case 0:
16102b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
16112b4ed282SShri Abhyankar         break;
16122b4ed282SShri Abhyankar       case 1:
16132b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
16142b4ed282SShri Abhyankar         break;
16152b4ed282SShri Abhyankar       case 2:
16162b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
16172b4ed282SShri Abhyankar         break;
16182b4ed282SShri Abhyankar       case 3:
16192b4ed282SShri Abhyankar         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
16202b4ed282SShri Abhyankar         break;
16212b4ed282SShri Abhyankar       }
16222b4ed282SShri Abhyankar     }
16232b4ed282SShri Abhyankar   ierr = PetscOptionsTail();CHKERRQ(ierr);
16242b4ed282SShri Abhyankar   PetscFunctionReturn(0);
16252b4ed282SShri Abhyankar }
16262b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
16272b4ed282SShri Abhyankar /*MC
16282b4ed282SShri Abhyankar       SNESVI - Semismooth newton method based nonlinear solver that uses a line search
16292b4ed282SShri Abhyankar 
16302b4ed282SShri Abhyankar    Options Database:
16312b4ed282SShri Abhyankar +   -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search
16322b4ed282SShri Abhyankar .   -snes_vi_alpha <alpha> - Sets alpha
16332b4ed282SShri 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)
16342b4ed282SShri Abhyankar .   -snes_vi_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
16352b4ed282SShri Abhyankar -   -snes_vi_monitor - print information about progress of line searches
16362b4ed282SShri Abhyankar 
16372b4ed282SShri Abhyankar 
16382b4ed282SShri Abhyankar    Level: beginner
16392b4ed282SShri Abhyankar 
16402b4ed282SShri Abhyankar .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
16412b4ed282SShri Abhyankar            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
16422b4ed282SShri Abhyankar            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
16432b4ed282SShri Abhyankar 
16442b4ed282SShri Abhyankar M*/
16452b4ed282SShri Abhyankar EXTERN_C_BEGIN
16462b4ed282SShri Abhyankar #undef __FUNCT__
16472b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI"
16482b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes)
16492b4ed282SShri Abhyankar {
16502b4ed282SShri Abhyankar   PetscErrorCode ierr;
16512b4ed282SShri Abhyankar   SNES_VI      *vi;
16522b4ed282SShri Abhyankar 
16532b4ed282SShri Abhyankar   PetscFunctionBegin;
16542b4ed282SShri Abhyankar   snes->ops->setup	     = SNESSetUp_VI;
165569c03620SShri Abhyankar   snes->ops->solve	     = SNESSolveVI_SS;
16562b4ed282SShri Abhyankar   snes->ops->destroy	     = SNESDestroy_VI;
16572b4ed282SShri Abhyankar   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
16582b4ed282SShri Abhyankar   snes->ops->view            = SNESView_VI;
16592b4ed282SShri Abhyankar   snes->ops->converged       = SNESDefaultConverged_VI;
16602b4ed282SShri Abhyankar 
16612b4ed282SShri Abhyankar   ierr                   = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
16622b4ed282SShri Abhyankar   snes->data    	 = (void*)vi;
16632b4ed282SShri Abhyankar   vi->alpha		 = 1.e-4;
16642b4ed282SShri Abhyankar   vi->maxstep		 = 1.e8;
16652b4ed282SShri Abhyankar   vi->minlambda         = 1.e-12;
16662b4ed282SShri Abhyankar   vi->LineSearch        = SNESLineSearchCubic_VI;
16672b4ed282SShri Abhyankar   vi->lsP               = PETSC_NULL;
16682b4ed282SShri Abhyankar   vi->postcheckstep     = PETSC_NULL;
16692b4ed282SShri Abhyankar   vi->postcheck         = PETSC_NULL;
16702b4ed282SShri Abhyankar   vi->precheckstep      = PETSC_NULL;
16712b4ed282SShri Abhyankar   vi->precheck          = PETSC_NULL;
16722b4ed282SShri Abhyankar   vi->const_tol         =  2.2204460492503131e-16;
16732b4ed282SShri Abhyankar   vi->computessfunction = ComputeFischerFunction;
16742b4ed282SShri Abhyankar 
16752b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
16762b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
16772b4ed282SShri Abhyankar 
16782b4ed282SShri Abhyankar   PetscFunctionReturn(0);
16792b4ed282SShri Abhyankar }
16802b4ed282SShri Abhyankar EXTERN_C_END
1681