xref: /petsc/src/snes/impls/vi/vi.c (revision 51defd01c63820b7cf54c8e98893c7df4640d693)
12b4ed282SShri Abhyankar #define PETSCSNES_DLL
22b4ed282SShri Abhyankar 
3d29d83b2SLisandro Dalcin #include "../src/snes/impls/vi/viimpl.h" /*I "petscsnes.h" I*/
4408da460SShri Abhyankar #include "../include/private/kspimpl.h"
52b4ed282SShri Abhyankar 
69308c137SBarry Smith #undef __FUNCT__
79308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
87087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
99308c137SBarry Smith {
109308c137SBarry Smith   PetscErrorCode          ierr;
119308c137SBarry Smith   SNES_VI                 *vi = (SNES_VI*)snes->data;
129308c137SBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;
139308c137SBarry Smith   const PetscScalar       *x,*xl,*xu,*f;
1409573ac7SBarry Smith   PetscInt                i,n,act = 0;
159308c137SBarry Smith   PetscReal               rnorm,fnorm;
169308c137SBarry Smith 
179308c137SBarry Smith   PetscFunctionBegin;
189308c137SBarry Smith   if (!dummy) {
199308c137SBarry Smith     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
209308c137SBarry Smith   }
219308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
229308c137SBarry Smith   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
239308c137SBarry Smith   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
249308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
253f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
269308c137SBarry Smith 
279308c137SBarry Smith   rnorm = 0.0;
289308c137SBarry Smith   for (i=0; i<n; i++) {
297b2c10efSShri Abhyankar     if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) rnorm += f[i]*f[i];
3009573ac7SBarry Smith     else act++;
319308c137SBarry Smith   }
323f731af5SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
339308c137SBarry Smith   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
349308c137SBarry Smith   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
359308c137SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
369308c137SBarry Smith   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
379308c137SBarry Smith   fnorm = sqrt(fnorm);
3809573ac7SBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES VI Function norm %14.12e Active constraints %D\n",its,fnorm,act);CHKERRQ(ierr);
399308c137SBarry Smith   if (!dummy) {
409308c137SBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
419308c137SBarry Smith   }
429308c137SBarry Smith   PetscFunctionReturn(0);
439308c137SBarry Smith }
449308c137SBarry Smith 
452b4ed282SShri Abhyankar /*
462b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
472b4ed282SShri 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
482b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
492b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
502b4ed282SShri Abhyankar */
512b4ed282SShri Abhyankar #undef __FUNCT__
522b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
53ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
542b4ed282SShri Abhyankar {
552b4ed282SShri Abhyankar   PetscReal      a1;
562b4ed282SShri Abhyankar   PetscErrorCode ierr;
57ace3abfcSBarry Smith   PetscBool     hastranspose;
582b4ed282SShri Abhyankar 
592b4ed282SShri Abhyankar   PetscFunctionBegin;
602b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
612b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
622b4ed282SShri Abhyankar   if (hastranspose) {
632b4ed282SShri Abhyankar     /* Compute || J^T F|| */
642b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
652b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
662b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
672b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
682b4ed282SShri Abhyankar   } else {
692b4ed282SShri Abhyankar     Vec         work;
702b4ed282SShri Abhyankar     PetscScalar result;
712b4ed282SShri Abhyankar     PetscReal   wnorm;
722b4ed282SShri Abhyankar 
732b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
742b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
752b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
762b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
772b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
782b4ed282SShri Abhyankar     ierr = VecDestroy(work);CHKERRQ(ierr);
792b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
802b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
812b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
822b4ed282SShri Abhyankar   }
832b4ed282SShri Abhyankar   PetscFunctionReturn(0);
842b4ed282SShri Abhyankar }
852b4ed282SShri Abhyankar 
862b4ed282SShri Abhyankar /*
872b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
882b4ed282SShri Abhyankar */
892b4ed282SShri Abhyankar #undef __FUNCT__
902b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
912b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
922b4ed282SShri Abhyankar {
932b4ed282SShri Abhyankar   PetscReal      a1,a2;
942b4ed282SShri Abhyankar   PetscErrorCode ierr;
95ace3abfcSBarry Smith   PetscBool     hastranspose;
962b4ed282SShri Abhyankar 
972b4ed282SShri Abhyankar   PetscFunctionBegin;
982b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
992b4ed282SShri Abhyankar   if (hastranspose) {
1002b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1012b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1022b4ed282SShri Abhyankar 
1032b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1042b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
1052b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
1062b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
1072b4ed282SShri Abhyankar     if (a1 != 0.0) {
1082b4ed282SShri Abhyankar       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
1092b4ed282SShri Abhyankar     }
1102b4ed282SShri Abhyankar   }
1112b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1122b4ed282SShri Abhyankar }
1132b4ed282SShri Abhyankar 
1142b4ed282SShri Abhyankar /*
1152b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
1162b4ed282SShri Abhyankar 
1172b4ed282SShri Abhyankar   Notes:
1182b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1192b4ed282SShri Abhyankar   merit < abstol
1202b4ed282SShri Abhyankar   merit < rtol*merit_initial
1212b4ed282SShri Abhyankar */
1222b4ed282SShri Abhyankar #undef __FUNCT__
1232b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
1247fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
1252b4ed282SShri Abhyankar {
1262b4ed282SShri Abhyankar   PetscErrorCode ierr;
1272b4ed282SShri Abhyankar 
1282b4ed282SShri Abhyankar   PetscFunctionBegin;
1292b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1302b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
1312b4ed282SShri Abhyankar 
1322b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1332b4ed282SShri Abhyankar 
1342b4ed282SShri Abhyankar   if (!it) {
1352b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1367fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
1372b4ed282SShri Abhyankar   }
1387fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
1392b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
1402b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1417fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
1427fe79bc4SShri Abhyankar     ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr);
1432b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
1442b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
1452b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
1462b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
1472b4ed282SShri Abhyankar   }
1482b4ed282SShri Abhyankar 
1492b4ed282SShri Abhyankar   if (it && !*reason) {
1507fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
1517fe79bc4SShri Abhyankar       ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr);
1522b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
1532b4ed282SShri Abhyankar     }
1542b4ed282SShri Abhyankar   }
1552b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1562b4ed282SShri Abhyankar }
1572b4ed282SShri Abhyankar 
1582b4ed282SShri Abhyankar /*
1592b4ed282SShri Abhyankar   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
1602b4ed282SShri Abhyankar 
1612b4ed282SShri Abhyankar   Input Parameter:
1622b4ed282SShri Abhyankar . phi - the semismooth function
1632b4ed282SShri Abhyankar 
1642b4ed282SShri Abhyankar   Output Parameter:
1652b4ed282SShri Abhyankar . merit - the merit function
1662b4ed282SShri Abhyankar . phinorm - ||phi||
1672b4ed282SShri Abhyankar 
1682b4ed282SShri Abhyankar   Notes:
1692b4ed282SShri Abhyankar   The merit function for the mixed complementarity problem is defined as
1702b4ed282SShri Abhyankar      merit = 0.5*phi^T*phi
1712b4ed282SShri Abhyankar */
1722b4ed282SShri Abhyankar #undef __FUNCT__
1732b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction"
1742b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
1752b4ed282SShri Abhyankar {
1762b4ed282SShri Abhyankar   PetscErrorCode ierr;
1772b4ed282SShri Abhyankar 
1782b4ed282SShri Abhyankar   PetscFunctionBegin;
1792b4ed282SShri Abhyankar   ierr = VecNormBegin(phi,NORM_2,phinorm);
1802b4ed282SShri Abhyankar   ierr = VecNormEnd(phi,NORM_2,phinorm);
1812b4ed282SShri Abhyankar 
1822b4ed282SShri Abhyankar   *merit = 0.5*(*phinorm)*(*phinorm);
1832b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1842b4ed282SShri Abhyankar }
1852b4ed282SShri Abhyankar 
1864c21c6cdSBarry Smith static inline PetscScalar Phi(PetscScalar a,PetscScalar b)
1874c21c6cdSBarry Smith {
1884c21c6cdSBarry Smith   return a + b - sqrt(a*a + b*b);
1894c21c6cdSBarry Smith }
1904c21c6cdSBarry Smith 
1914c21c6cdSBarry Smith static inline PetscScalar DPhi(PetscScalar a,PetscScalar b)
1924c21c6cdSBarry Smith {
1935559b345SBarry Smith   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return  1.0 - a/ sqrt(a*a + b*b);
1945559b345SBarry Smith   else return .5;
1954c21c6cdSBarry Smith }
1964c21c6cdSBarry Smith 
1972b4ed282SShri Abhyankar /*
1981f79c099SShri Abhyankar    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
1992b4ed282SShri Abhyankar 
2002b4ed282SShri Abhyankar    Input Parameters:
2012b4ed282SShri Abhyankar .  snes - the SNES context
2022b4ed282SShri Abhyankar .  x - current iterate
2032b4ed282SShri Abhyankar .  functx - user defined function context
2042b4ed282SShri Abhyankar 
2052b4ed282SShri Abhyankar    Output Parameters:
2062b4ed282SShri Abhyankar .  phi - Semismooth function
2072b4ed282SShri Abhyankar 
2082b4ed282SShri Abhyankar */
2092b4ed282SShri Abhyankar #undef __FUNCT__
2101f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction"
2111f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
2122b4ed282SShri Abhyankar {
2132b4ed282SShri Abhyankar   PetscErrorCode  ierr;
2142b4ed282SShri Abhyankar   SNES_VI         *vi = (SNES_VI*)snes->data;
2152b4ed282SShri Abhyankar   Vec             Xl = vi->xl,Xu = vi->xu,F = snes->vec_func;
2164c21c6cdSBarry Smith   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u;
2172b4ed282SShri Abhyankar   PetscInt        i,nlocal;
2182b4ed282SShri Abhyankar 
2192b4ed282SShri Abhyankar   PetscFunctionBegin;
2202b4ed282SShri Abhyankar   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
2212b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2222b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
2232b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
2242b4ed282SShri Abhyankar   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
2252b4ed282SShri Abhyankar   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
2262b4ed282SShri Abhyankar   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
2272b4ed282SShri Abhyankar 
2282b4ed282SShri Abhyankar   for (i=0;i < nlocal;i++) {
2295559b345SBarry Smith     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { /* no constraints on variable */
2304c21c6cdSBarry Smith       phi_arr[i] = f_arr[i];
2315559b345SBarry Smith     } else if (l[i] <= PETSC_VI_NINF) {                      /* upper bound on variable only */
2324c21c6cdSBarry Smith       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
2335559b345SBarry Smith     } else if (u[i] >= PETSC_VI_INF) {                       /* lower bound on variable only */
2344c21c6cdSBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
2355559b345SBarry Smith     } else if (l[i] == u[i]) {
2362b4ed282SShri Abhyankar       phi_arr[i] = l[i] - x_arr[i];
2375559b345SBarry Smith     } else {                                                /* both bounds on variable */
2384c21c6cdSBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
2392b4ed282SShri Abhyankar     }
2402b4ed282SShri Abhyankar   }
2412b4ed282SShri Abhyankar 
2422b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
2432b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
2442b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
2452b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
2462b4ed282SShri Abhyankar   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
2472b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2482b4ed282SShri Abhyankar }
2492b4ed282SShri Abhyankar 
2502b4ed282SShri Abhyankar /*
251a79edbefSShri Abhyankar    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
252a79edbefSShri Abhyankar                                           the semismooth jacobian.
2532b4ed282SShri Abhyankar */
2542b4ed282SShri Abhyankar #undef __FUNCT__
255a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
256a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
2572b4ed282SShri Abhyankar {
2582b4ed282SShri Abhyankar   PetscErrorCode ierr;
2592b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*)snes->data;
2604c21c6cdSBarry Smith   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
2612b4ed282SShri Abhyankar   PetscInt       i,nlocal;
2622b4ed282SShri Abhyankar 
2632b4ed282SShri Abhyankar   PetscFunctionBegin;
2642b4ed282SShri Abhyankar 
2652b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
2662b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
2672b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr);
2682b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr);
269a79edbefSShri Abhyankar   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
270a79edbefSShri Abhyankar   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
2712b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2724c21c6cdSBarry Smith 
2732b4ed282SShri Abhyankar   for (i=0;i< nlocal;i++) {
2745559b345SBarry Smith     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {/* no constraints on variable */
2754c21c6cdSBarry Smith       da[i] = 0;
2762b4ed282SShri Abhyankar       db[i] = 1;
2775559b345SBarry Smith     } else if (l[i] <= PETSC_VI_NINF) {                     /* upper bound on variable only */
2784c21c6cdSBarry Smith       da[i] = DPhi(u[i] - x[i], -f[i]);
2794c21c6cdSBarry Smith       db[i] = DPhi(-f[i],u[i] - x[i]);
2805559b345SBarry Smith     } else if (u[i] >= PETSC_VI_INF) {                      /* lower bound on variable only */
2815559b345SBarry Smith       da[i] = DPhi(x[i] - l[i], f[i]);
2825559b345SBarry Smith       db[i] = DPhi(f[i],x[i] - l[i]);
2835559b345SBarry Smith     } else if (l[i] == u[i]) {                              /* fixed variable */
2844c21c6cdSBarry Smith       da[i] = 1;
2852b4ed282SShri Abhyankar       db[i] = 0;
2865559b345SBarry Smith     } else {                                                /* upper and lower bounds on variable */
2874c21c6cdSBarry Smith       da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
2884c21c6cdSBarry Smith       db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
2894c21c6cdSBarry Smith       da2 = DPhi(u[i] - x[i], -f[i]);
2904c21c6cdSBarry Smith       db2 = DPhi(-f[i],u[i] - x[i]);
2914c21c6cdSBarry Smith       da[i] = da1 + db1*da2;
2924c21c6cdSBarry Smith       db[i] = db1*db2;
2932b4ed282SShri Abhyankar     }
2942b4ed282SShri Abhyankar   }
2955559b345SBarry Smith 
2962b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
2972b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
2982b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr);
2992b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr);
300a79edbefSShri Abhyankar   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
301a79edbefSShri Abhyankar   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
302a79edbefSShri Abhyankar   PetscFunctionReturn(0);
303a79edbefSShri Abhyankar }
304a79edbefSShri Abhyankar 
305a79edbefSShri Abhyankar /*
306a79edbefSShri 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.
307a79edbefSShri Abhyankar 
308a79edbefSShri Abhyankar    Input Parameters:
309a79edbefSShri Abhyankar .  Da       - Diagonal shift vector for the semismooth jacobian.
310a79edbefSShri Abhyankar .  Db       - Row scaling vector for the semismooth jacobian.
311a79edbefSShri Abhyankar 
312a79edbefSShri Abhyankar    Output Parameters:
313a79edbefSShri Abhyankar .  jac      - semismooth jacobian
314a79edbefSShri Abhyankar .  jac_pre  - optional preconditioning matrix
315a79edbefSShri Abhyankar 
316a79edbefSShri Abhyankar    Notes:
317a79edbefSShri Abhyankar    The semismooth jacobian matrix is given by
318a79edbefSShri Abhyankar    jac = Da + Db*jacfun
319a79edbefSShri Abhyankar    where Db is the row scaling matrix stored as a vector,
320a79edbefSShri Abhyankar          Da is the diagonal perturbation matrix stored as a vector
321a79edbefSShri Abhyankar    and   jacfun is the jacobian of the original nonlinear function.
322a79edbefSShri Abhyankar */
323a79edbefSShri Abhyankar #undef __FUNCT__
324a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian"
325a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
326a79edbefSShri Abhyankar {
327a79edbefSShri Abhyankar   PetscErrorCode ierr;
328a79edbefSShri Abhyankar 
3292b4ed282SShri Abhyankar   /* Do row scaling  and add diagonal perturbation */
330a79edbefSShri Abhyankar   ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr);
331a79edbefSShri Abhyankar   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
332a79edbefSShri Abhyankar   if (jac != jac_pre) { /* If jac and jac_pre are different */
333a79edbefSShri Abhyankar     ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL);
334a79edbefSShri Abhyankar     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
3352b4ed282SShri Abhyankar   }
3362b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3372b4ed282SShri Abhyankar }
3382b4ed282SShri Abhyankar 
3392b4ed282SShri Abhyankar /*
340ee657d29SShri Abhyankar    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
341ee657d29SShri Abhyankar 
342ee657d29SShri Abhyankar    Input Parameters:
343ee657d29SShri Abhyankar    phi - semismooth function.
344ee657d29SShri Abhyankar    H   - semismooth jacobian
345ee657d29SShri Abhyankar 
346ee657d29SShri Abhyankar    Output Parameters:
347ee657d29SShri Abhyankar    dpsi - merit function gradient
348ee657d29SShri Abhyankar 
349ee657d29SShri Abhyankar    Notes:
350ee657d29SShri Abhyankar   The merit function gradient is computed as follows
351ee657d29SShri Abhyankar         dpsi = H^T*phi
352ee657d29SShri Abhyankar */
353ee657d29SShri Abhyankar #undef __FUNCT__
354ee657d29SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
355ee657d29SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
356ee657d29SShri Abhyankar {
357ee657d29SShri Abhyankar   PetscErrorCode ierr;
358ee657d29SShri Abhyankar 
359ee657d29SShri Abhyankar   PetscFunctionBegin;
360ee657d29SShri Abhyankar   ierr = MatMultTranspose(H,phi,dpsi);
361ee657d29SShri Abhyankar   PetscFunctionReturn(0);
362ee657d29SShri Abhyankar }
363ee657d29SShri Abhyankar 
364c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
365c1a5e521SShri Abhyankar /*
366c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
367c1a5e521SShri Abhyankar 
368c1a5e521SShri Abhyankar    Input Parameters:
369c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
370c1a5e521SShri Abhyankar 
371c1a5e521SShri Abhyankar    Output Parameters:
372c1a5e521SShri Abhyankar .  X - Bound projected X
373c1a5e521SShri Abhyankar 
374c1a5e521SShri Abhyankar */
375c1a5e521SShri Abhyankar 
376c1a5e521SShri Abhyankar #undef __FUNCT__
377c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
378c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
379c1a5e521SShri Abhyankar {
380c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
381c1a5e521SShri Abhyankar   SNES_VI           *vi = (SNES_VI*)snes->data;
382c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
383c1a5e521SShri Abhyankar   PetscScalar       *x;
384c1a5e521SShri Abhyankar   PetscInt          i,n;
385c1a5e521SShri Abhyankar 
386c1a5e521SShri Abhyankar   PetscFunctionBegin;
387c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
388c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
389c1a5e521SShri Abhyankar   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
390c1a5e521SShri Abhyankar   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
391c1a5e521SShri Abhyankar 
392c1a5e521SShri Abhyankar   for(i = 0;i<n;i++) {
393c1a5e521SShri Abhyankar     if (x[i] < xl[i]) x[i] = xl[i];
394c1a5e521SShri Abhyankar     else if (x[i] > xu[i]) x[i] = xu[i];
395c1a5e521SShri Abhyankar   }
396c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
397c1a5e521SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
398c1a5e521SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
399c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
400c1a5e521SShri Abhyankar }
401c1a5e521SShri Abhyankar 
4022b4ed282SShri Abhyankar /*  --------------------------------------------------------------------
4032b4ed282SShri Abhyankar 
4042b4ed282SShri Abhyankar      This file implements a semismooth truncated Newton method with a line search,
4052b4ed282SShri Abhyankar      for solving a system of nonlinear equations in complementarity form, using the KSP, Vec,
4062b4ed282SShri Abhyankar      and Mat interfaces for linear solvers, vectors, and matrices,
4072b4ed282SShri Abhyankar      respectively.
4082b4ed282SShri Abhyankar 
4092b4ed282SShri Abhyankar      The following basic routines are required for each nonlinear solver:
4102b4ed282SShri Abhyankar           SNESCreate_XXX()          - Creates a nonlinear solver context
4112b4ed282SShri Abhyankar           SNESSetFromOptions_XXX()  - Sets runtime options
4122b4ed282SShri Abhyankar           SNESSolve_XXX()           - Solves the nonlinear system
4132b4ed282SShri Abhyankar           SNESDestroy_XXX()         - Destroys the nonlinear solver context
4142b4ed282SShri Abhyankar      The suffix "_XXX" denotes a particular implementation, in this case
4152b4ed282SShri Abhyankar      we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving
4162b4ed282SShri Abhyankar      systems of nonlinear equations with a line search (LS) method.
4172b4ed282SShri Abhyankar      These routines are actually called via the common user interface
4182b4ed282SShri Abhyankar      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
4192b4ed282SShri Abhyankar      SNESDestroy(), so the application code interface remains identical
4202b4ed282SShri Abhyankar      for all nonlinear solvers.
4212b4ed282SShri Abhyankar 
4222b4ed282SShri Abhyankar      Another key routine is:
4232b4ed282SShri Abhyankar           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
4242b4ed282SShri Abhyankar      by setting data structures and options.   The interface routine SNESSetUp()
4252b4ed282SShri Abhyankar      is not usually called directly by the user, but instead is called by
4262b4ed282SShri Abhyankar      SNESSolve() if necessary.
4272b4ed282SShri Abhyankar 
4282b4ed282SShri Abhyankar      Additional basic routines are:
4292b4ed282SShri Abhyankar           SNESView_XXX()            - Prints details of runtime options that
4302b4ed282SShri Abhyankar                                       have actually been used.
4312b4ed282SShri Abhyankar      These are called by application codes via the interface routines
4322b4ed282SShri Abhyankar      SNESView().
4332b4ed282SShri Abhyankar 
4342b4ed282SShri Abhyankar      The various types of solvers (preconditioners, Krylov subspace methods,
4352b4ed282SShri Abhyankar      nonlinear solvers, timesteppers) are all organized similarly, so the
4362b4ed282SShri Abhyankar      above description applies to these categories also.
4372b4ed282SShri Abhyankar 
4382b4ed282SShri Abhyankar     -------------------------------------------------------------------- */
4392b4ed282SShri Abhyankar /*
44069c03620SShri Abhyankar    SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton
4412b4ed282SShri Abhyankar    method using a line search.
4422b4ed282SShri Abhyankar 
4432b4ed282SShri Abhyankar    Input Parameters:
4442b4ed282SShri Abhyankar .  snes - the SNES context
4452b4ed282SShri Abhyankar 
4462b4ed282SShri Abhyankar    Output Parameter:
4472b4ed282SShri Abhyankar .  outits - number of iterations until termination
4482b4ed282SShri Abhyankar 
4492b4ed282SShri Abhyankar    Application Interface Routine: SNESSolve()
4502b4ed282SShri Abhyankar 
4512b4ed282SShri Abhyankar    Notes:
4522b4ed282SShri Abhyankar    This implements essentially a semismooth Newton method with a
453*51defd01SShri Abhyankar    line search. The default line search does not do any line seach
454*51defd01SShri Abhyankar    but rather takes a full newton step.
4552b4ed282SShri Abhyankar */
4562b4ed282SShri Abhyankar #undef __FUNCT__
45769c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS"
45869c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes)
4592b4ed282SShri Abhyankar {
4602b4ed282SShri Abhyankar   SNES_VI            *vi = (SNES_VI*)snes->data;
4612b4ed282SShri Abhyankar   PetscErrorCode     ierr;
4622b4ed282SShri Abhyankar   PetscInt           maxits,i,lits;
4633b336b49SShri Abhyankar   PetscBool          lssucceed;
4642b4ed282SShri Abhyankar   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
4652b4ed282SShri Abhyankar   PetscReal          gnorm,xnorm=0,ynorm;
4662b4ed282SShri Abhyankar   Vec                Y,X,F,G,W;
4672b4ed282SShri Abhyankar   KSPConvergedReason kspreason;
4682b4ed282SShri Abhyankar 
4692b4ed282SShri Abhyankar   PetscFunctionBegin;
4709ce7406fSBarry Smith   vi->computeuserfunction    = snes->ops->computefunction;
4719ce7406fSBarry Smith   snes->ops->computefunction = SNESVIComputeFunction;
4729ce7406fSBarry Smith 
4732b4ed282SShri Abhyankar   snes->numFailures            = 0;
4742b4ed282SShri Abhyankar   snes->numLinearSolveFailures = 0;
4752b4ed282SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
4762b4ed282SShri Abhyankar 
4772b4ed282SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
4782b4ed282SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
4792b4ed282SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
4802b4ed282SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
4812b4ed282SShri Abhyankar   G		= snes->work[1];
4822b4ed282SShri Abhyankar   W		= snes->work[2];
4832b4ed282SShri Abhyankar 
4842b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
4852b4ed282SShri Abhyankar   snes->iter = 0;
4862b4ed282SShri Abhyankar   snes->norm = 0.0;
4872b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
4882b4ed282SShri Abhyankar 
4897fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
4902b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
4912b4ed282SShri Abhyankar   if (snes->domainerror) {
4922b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
4939ce7406fSBarry Smith     snes->ops->computefunction = vi->computeuserfunction;
4942b4ed282SShri Abhyankar     PetscFunctionReturn(0);
4952b4ed282SShri Abhyankar   }
4962b4ed282SShri Abhyankar    /* Compute Merit function */
4972b4ed282SShri Abhyankar   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
4982b4ed282SShri Abhyankar 
4992b4ed282SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
5002b4ed282SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
5012b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5022b4ed282SShri Abhyankar 
5032b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5042b4ed282SShri Abhyankar   snes->norm = vi->phinorm;
5052b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5062b4ed282SShri Abhyankar   SNESLogConvHistory(snes,vi->phinorm,0);
5072b4ed282SShri Abhyankar   SNESMonitor(snes,0,vi->phinorm);
5082b4ed282SShri Abhyankar 
5092b4ed282SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
5102b4ed282SShri Abhyankar   snes->ttol = vi->phinorm*snes->rtol;
5112b4ed282SShri Abhyankar   /* test convergence */
5122b4ed282SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
5139ce7406fSBarry Smith   if (snes->reason) {
5149ce7406fSBarry Smith     snes->ops->computefunction = vi->computeuserfunction;
5159ce7406fSBarry Smith     PetscFunctionReturn(0);
5169ce7406fSBarry Smith   }
5172b4ed282SShri Abhyankar 
5182b4ed282SShri Abhyankar   for (i=0; i<maxits; i++) {
5192b4ed282SShri Abhyankar 
5202b4ed282SShri Abhyankar     /* Call general purpose update function */
5212b4ed282SShri Abhyankar     if (snes->ops->update) {
5222b4ed282SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
5232b4ed282SShri Abhyankar     }
5242b4ed282SShri Abhyankar 
5252b4ed282SShri Abhyankar     /* Solve J Y = Phi, where J is the semismooth jacobian */
526a79edbefSShri Abhyankar     /* Get the nonlinear function jacobian */
5272b4ed282SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
528a79edbefSShri Abhyankar     /* Get the diagonal shift and row scaling vectors */
529a79edbefSShri Abhyankar     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
530a79edbefSShri Abhyankar     /* Compute the semismooth jacobian */
531a79edbefSShri Abhyankar     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
532ee657d29SShri Abhyankar     /* Compute the merit function gradient */
533ee657d29SShri Abhyankar     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
5342b4ed282SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
5352b4ed282SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
5362b4ed282SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
5373b336b49SShri Abhyankar 
5383b336b49SShri Abhyankar     if (kspreason < 0) {
5392b4ed282SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
5402b4ed282SShri 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);
5412b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
5422b4ed282SShri Abhyankar         break;
5432b4ed282SShri Abhyankar       }
5442b4ed282SShri Abhyankar     }
5452b4ed282SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
5462b4ed282SShri Abhyankar     snes->linear_its += lits;
5472b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
5482b4ed282SShri Abhyankar     /*
5492b4ed282SShri Abhyankar     if (vi->precheckstep) {
550ace3abfcSBarry Smith       PetscBool changed_y = PETSC_FALSE;
5512b4ed282SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
5522b4ed282SShri Abhyankar     }
5532b4ed282SShri Abhyankar 
5542b4ed282SShri Abhyankar     if (PetscLogPrintInfo){
5552b4ed282SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
5562b4ed282SShri Abhyankar     }
5572b4ed282SShri Abhyankar     */
5582b4ed282SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
5592b4ed282SShri Abhyankar          Y <- X - lambda*Y
5602b4ed282SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
5612b4ed282SShri Abhyankar     */
5622b4ed282SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
5632b4ed282SShri Abhyankar     ynorm = 1; gnorm = vi->phinorm;
5642b4ed282SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
5652b4ed282SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
5662b4ed282SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
5672b4ed282SShri Abhyankar     if (snes->domainerror) {
5682b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
5699ce7406fSBarry Smith       snes->ops->computefunction = vi->computeuserfunction;
5702b4ed282SShri Abhyankar       PetscFunctionReturn(0);
5712b4ed282SShri Abhyankar     }
5722b4ed282SShri Abhyankar     if (!lssucceed) {
5732b4ed282SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
574ace3abfcSBarry Smith 	PetscBool ismin;
5752b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
5762b4ed282SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
5772b4ed282SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
5782b4ed282SShri Abhyankar         break;
5792b4ed282SShri Abhyankar       }
5802b4ed282SShri Abhyankar     }
5812b4ed282SShri Abhyankar     /* Update function and solution vectors */
5822b4ed282SShri Abhyankar     vi->phinorm = gnorm;
5832b4ed282SShri Abhyankar     vi->merit = 0.5*vi->phinorm*vi->phinorm;
5842b4ed282SShri Abhyankar     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
5852b4ed282SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
5862b4ed282SShri Abhyankar     /* Monitor convergence */
5872b4ed282SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5882b4ed282SShri Abhyankar     snes->iter = i+1;
5892b4ed282SShri Abhyankar     snes->norm = vi->phinorm;
5902b4ed282SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5912b4ed282SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
5922b4ed282SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
5932b4ed282SShri Abhyankar     /* Test for convergence, xnorm = || X || */
5942b4ed282SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
5952b4ed282SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
5962b4ed282SShri Abhyankar     if (snes->reason) break;
5972b4ed282SShri Abhyankar   }
5982b4ed282SShri Abhyankar   if (i == maxits) {
5992b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
6002b4ed282SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
6012b4ed282SShri Abhyankar   }
6029ce7406fSBarry Smith   snes->ops->computefunction = vi->computeuserfunction;
6032b4ed282SShri Abhyankar   PetscFunctionReturn(0);
6042b4ed282SShri Abhyankar }
60569c03620SShri Abhyankar 
60669c03620SShri Abhyankar #undef __FUNCT__
607d950fb63SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
608d950fb63SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec Xl, Vec Xu,IS* ISact,IS* ISinact)
609d950fb63SShri Abhyankar {
610d950fb63SShri Abhyankar   PetscErrorCode     ierr;
611d950fb63SShri Abhyankar   PetscInt           i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0;
612d950fb63SShri Abhyankar   PetscInt          *idx_act,*idx_inact,i1=0,i2=0;
613fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
614d950fb63SShri Abhyankar   Vec                F = snes->vec_func;
615d950fb63SShri Abhyankar 
616d950fb63SShri Abhyankar   PetscFunctionBegin;
617d950fb63SShri Abhyankar 
618d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
619d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
620fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
621fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
622fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
623fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
624d950fb63SShri Abhyankar   /* Compute the sizes of the active and inactive sets */
625d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
626fe835674SShri Abhyankar     if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) nloc_isinact++;
627fe835674SShri Abhyankar     else nloc_isact++;
628d950fb63SShri Abhyankar   }
629d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
630d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr);
631d950fb63SShri Abhyankar 
632d950fb63SShri Abhyankar   /* Creating the indexing arrays */
633d950fb63SShri Abhyankar   for(i=0; i < nlocal; i++) {
634fe835674SShri Abhyankar     if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) idx_inact[i2++] = ilow+i;
635fe835674SShri Abhyankar     else idx_act[i1++] = ilow+i;
636d950fb63SShri Abhyankar   }
637d950fb63SShri Abhyankar 
638d950fb63SShri Abhyankar   /* Create the index sets */
639d950fb63SShri Abhyankar   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr);
640d950fb63SShri Abhyankar   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr);
641d950fb63SShri Abhyankar 
642fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
643fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
644fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
645fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
646d950fb63SShri Abhyankar   ierr = PetscFree(idx_act);CHKERRQ(ierr);
647d950fb63SShri Abhyankar   ierr = PetscFree(idx_inact);CHKERRQ(ierr);
648d950fb63SShri Abhyankar   PetscFunctionReturn(0);
649d950fb63SShri Abhyankar }
650d950fb63SShri Abhyankar 
651dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
652dbd914b8SShri Abhyankar #undef __FUNCT__
653fe835674SShri Abhyankar #define __FUNCT__ "SNESVICreateSubVectors"
654fe835674SShri Abhyankar PetscErrorCode SNESVICreateSubVectors(SNES snes,PetscInt n,Vec* newv)
655dbd914b8SShri Abhyankar {
656dbd914b8SShri Abhyankar   PetscErrorCode ierr;
657dbd914b8SShri Abhyankar   Vec            v;
658dbd914b8SShri Abhyankar 
659dbd914b8SShri Abhyankar   PetscFunctionBegin;
660dbd914b8SShri Abhyankar   ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr);
661dbd914b8SShri Abhyankar   ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr);
662dbd914b8SShri Abhyankar   ierr = VecSetFromOptions(v);CHKERRQ(ierr);
663dbd914b8SShri Abhyankar   *newv = v;
664dbd914b8SShri Abhyankar 
665dbd914b8SShri Abhyankar   PetscFunctionReturn(0);
666dbd914b8SShri Abhyankar }
667dbd914b8SShri Abhyankar 
668732bb160SShri Abhyankar /* Resets the snes PC and KSP when the active set sizes change */
669732bb160SShri Abhyankar #undef __FUNCT__
670732bb160SShri Abhyankar #define __FUNCT__ "SNESVIResetPCandKSP"
671732bb160SShri Abhyankar PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
672732bb160SShri Abhyankar {
673732bb160SShri Abhyankar   PetscErrorCode         ierr;
674732bb160SShri Abhyankar   KSP                    kspnew,snesksp;
675732bb160SShri Abhyankar   PC                     pcnew;
676732bb160SShri Abhyankar   const MatSolverPackage stype;
677dbd914b8SShri Abhyankar 
678732bb160SShri Abhyankar   PetscFunctionBegin;
679732bb160SShri Abhyankar   /* The active and inactive set sizes have changed so need to create a new snes->ksp object */
680732bb160SShri Abhyankar   ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
681732bb160SShri Abhyankar   ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr);
682732bb160SShri Abhyankar   /* Copy over snes->ksp info */
683732bb160SShri Abhyankar   kspnew->pc_side = snesksp->pc_side;
684732bb160SShri Abhyankar   kspnew->rtol    = snesksp->rtol;
685732bb160SShri Abhyankar   kspnew->abstol    = snesksp->abstol;
686732bb160SShri Abhyankar   kspnew->max_it  = snesksp->max_it;
687732bb160SShri Abhyankar   ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr);
688732bb160SShri Abhyankar   ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr);
689732bb160SShri Abhyankar   ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr);
690732bb160SShri Abhyankar   ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
691732bb160SShri Abhyankar   ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr);
692732bb160SShri Abhyankar   ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr);
693732bb160SShri Abhyankar   ierr = KSPDestroy(snesksp);CHKERRQ(ierr);
694732bb160SShri Abhyankar   snes->ksp = kspnew;
695732bb160SShri Abhyankar   ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr);
696732bb160SShri Abhyankar   ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);
697732bb160SShri Abhyankar   PetscFunctionReturn(0);
698732bb160SShri Abhyankar }
69969c03620SShri Abhyankar 
70069c03620SShri Abhyankar 
701fe835674SShri Abhyankar #undef __FUNCT__
702fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
703fe835674SShri Abhyankar PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscScalar *fnorm)
704fe835674SShri Abhyankar {
705fe835674SShri Abhyankar   PetscErrorCode    ierr;
706fe835674SShri Abhyankar   SNES_VI           *vi = (SNES_VI*)snes->data;
707fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
708fe835674SShri Abhyankar   PetscInt          i,n;
709fe835674SShri Abhyankar   PetscReal         rnorm;
710fe835674SShri Abhyankar 
711fe835674SShri Abhyankar   PetscFunctionBegin;
712fe835674SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
713fe835674SShri Abhyankar   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
714fe835674SShri Abhyankar   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
715fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
716fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
717fe835674SShri Abhyankar   rnorm = 0.0;
718fe835674SShri Abhyankar   for (i=0; i<n; i++) {
719fe835674SShri Abhyankar     if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) rnorm += f[i]*f[i];
720fe835674SShri Abhyankar   }
721fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
722fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
723fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
724fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
725fe835674SShri Abhyankar   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
726fe835674SShri Abhyankar   *fnorm = sqrt(*fnorm);
727fe835674SShri Abhyankar   PetscFunctionReturn(0);
728fe835674SShri Abhyankar }
729fe835674SShri Abhyankar 
730fe835674SShri Abhyankar /* Variational Inequality solver using reduce space method. No semismooth algorithm is
731fe835674SShri Abhyankar    implemented in this algorithm. It basically identifies the active variables and does
732fe835674SShri Abhyankar    a linear solve on the inactive variables. */
733d950fb63SShri Abhyankar #undef __FUNCT__
734d950fb63SShri Abhyankar #define __FUNCT__ "SNESSolveVI_RS"
735d950fb63SShri Abhyankar PetscErrorCode SNESSolveVI_RS(SNES snes)
736d950fb63SShri Abhyankar {
737d950fb63SShri Abhyankar   SNES_VI          *vi = (SNES_VI*)snes->data;
738d950fb63SShri Abhyankar   PetscErrorCode    ierr;
73964221577SShri Abhyankar   PetscInt          maxits,i,lits,Nis_inact;
740d950fb63SShri Abhyankar   PetscBool         lssucceed;
741d950fb63SShri Abhyankar   MatStructure      flg = DIFFERENT_NONZERO_PATTERN;
742fe835674SShri Abhyankar   PetscReal         fnorm,gnorm,xnorm=0,ynorm;
743d950fb63SShri Abhyankar   Vec                Y,X,F,G,W;
744d950fb63SShri Abhyankar   KSPConvergedReason kspreason;
745d950fb63SShri Abhyankar 
746d950fb63SShri Abhyankar   PetscFunctionBegin;
747d950fb63SShri Abhyankar   snes->numFailures            = 0;
748d950fb63SShri Abhyankar   snes->numLinearSolveFailures = 0;
749d950fb63SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
750d950fb63SShri Abhyankar 
751d950fb63SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
752d950fb63SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
753d950fb63SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
754d950fb63SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
755d950fb63SShri Abhyankar   G		= snes->work[1];
756d950fb63SShri Abhyankar   W		= snes->work[2];
757d950fb63SShri Abhyankar 
75864221577SShri Abhyankar   Nis_inact = F->map->N;
759d950fb63SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
760d950fb63SShri Abhyankar   snes->iter = 0;
761d950fb63SShri Abhyankar   snes->norm = 0.0;
762d950fb63SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
763d950fb63SShri Abhyankar 
7647fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
765fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
766d950fb63SShri Abhyankar   if (snes->domainerror) {
767d950fb63SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
768d950fb63SShri Abhyankar     PetscFunctionReturn(0);
769d950fb63SShri Abhyankar   }
770fe835674SShri Abhyankar   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
771d950fb63SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
772d950fb63SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
773fe835674SShri Abhyankar   if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
774d950fb63SShri Abhyankar 
775d950fb63SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
776fe835674SShri Abhyankar   snes->norm = fnorm;
777d950fb63SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
778fe835674SShri Abhyankar   SNESLogConvHistory(snes,fnorm,0);
779fe835674SShri Abhyankar   SNESMonitor(snes,0,fnorm);
780d950fb63SShri Abhyankar 
781d950fb63SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
782fe835674SShri Abhyankar   snes->ttol = fnorm*snes->rtol;
783d950fb63SShri Abhyankar   /* test convergence */
784fe835674SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
785d950fb63SShri Abhyankar   if (snes->reason) PetscFunctionReturn(0);
786d950fb63SShri Abhyankar 
787d950fb63SShri Abhyankar   for (i=0; i<maxits; i++) {
788d950fb63SShri Abhyankar 
789d950fb63SShri Abhyankar     IS                 IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
790d950fb63SShri Abhyankar     VecScatter         scat_act,scat_inact;
79164221577SShri Abhyankar     PetscInt           nis_act,nis_inact,Nis_inact_prev;
792fe835674SShri Abhyankar     Vec                Y_act,Y_inact,F_inact;
793d950fb63SShri Abhyankar     Mat                jac_inact_inact,prejac_inact_inact;
794d950fb63SShri Abhyankar 
795d950fb63SShri Abhyankar     /* Call general purpose update function */
796d950fb63SShri Abhyankar     if (snes->ops->update) {
797d950fb63SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
798d950fb63SShri Abhyankar     }
799d950fb63SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
800d950fb63SShri Abhyankar     /* Create active and inactive index sets */
801d950fb63SShri Abhyankar     ierr = SNESVICreateIndexSets_RS(snes,X,vi->xl,vi->xu,&IS_act,&IS_inact);CHKERRQ(ierr);
802d950fb63SShri Abhyankar 
80364221577SShri Abhyankar     Nis_inact_prev = Nis_inact;
804d950fb63SShri Abhyankar     /* Get sizes of active and inactive sets */
805d950fb63SShri Abhyankar     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
806d950fb63SShri Abhyankar     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
80764221577SShri Abhyankar     ierr = ISGetSize(IS_inact,&Nis_inact);CHKERRQ(ierr);
808d950fb63SShri Abhyankar 
809d950fb63SShri Abhyankar     ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr);
810d950fb63SShri Abhyankar 
811d950fb63SShri Abhyankar     /* Create active and inactive set vectors */
812fe835674SShri Abhyankar     ierr = SNESVICreateSubVectors(snes,nis_inact,&F_inact);CHKERRQ(ierr);
813fe835674SShri Abhyankar     ierr = SNESVICreateSubVectors(snes,nis_act,&Y_act);CHKERRQ(ierr);
814fe835674SShri Abhyankar     ierr = SNESVICreateSubVectors(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
815d950fb63SShri Abhyankar 
816d950fb63SShri Abhyankar     /* Create inactive set submatrices */
817d950fb63SShri Abhyankar     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
818d950fb63SShri Abhyankar 
819d950fb63SShri Abhyankar     /* Create scatter contexts */
820d950fb63SShri Abhyankar     ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
821d950fb63SShri Abhyankar     ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
822d950fb63SShri Abhyankar 
823d950fb63SShri Abhyankar     /* Do a vec scatter to active and inactive set vectors */
824fe835674SShri Abhyankar     ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
825fe835674SShri Abhyankar     ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
826d950fb63SShri Abhyankar 
827d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
828d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
829d950fb63SShri Abhyankar 
830d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
831d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
832d950fb63SShri Abhyankar 
833d950fb63SShri Abhyankar     /* Active set direction = 0 */
834d950fb63SShri Abhyankar     ierr = VecSet(Y_act,0);CHKERRQ(ierr);
835d950fb63SShri Abhyankar     if (snes->jacobian != snes->jacobian_pre) {
836d950fb63SShri Abhyankar       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
837d950fb63SShri Abhyankar     } else prejac_inact_inact = jac_inact_inact;
838d950fb63SShri Abhyankar 
83964221577SShri Abhyankar     if (Nis_inact != Nis_inact_prev) {
840732bb160SShri Abhyankar       ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr);
841d950fb63SShri Abhyankar     }
842d950fb63SShri Abhyankar 
843d950fb63SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
844fe835674SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr);
845d950fb63SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
846d950fb63SShri Abhyankar     if (kspreason < 0) {
847d950fb63SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
848d950fb63SShri 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);
849d950fb63SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
850d950fb63SShri Abhyankar         break;
851d950fb63SShri Abhyankar       }
852d950fb63SShri Abhyankar      }
853d950fb63SShri Abhyankar 
854d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
855d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
856d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
857d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
858d950fb63SShri Abhyankar 
859fe835674SShri Abhyankar     ierr = VecDestroy(F_inact);CHKERRQ(ierr);
860d950fb63SShri Abhyankar     ierr = VecDestroy(Y_act);CHKERRQ(ierr);
861d950fb63SShri Abhyankar     ierr = VecDestroy(Y_inact);CHKERRQ(ierr);
862d950fb63SShri Abhyankar     ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr);
863d950fb63SShri Abhyankar     ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr);
864d950fb63SShri Abhyankar     ierr = ISDestroy(IS_act);CHKERRQ(ierr);
865d950fb63SShri Abhyankar     ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
866d950fb63SShri Abhyankar     ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
867d950fb63SShri Abhyankar     if (snes->jacobian != snes->jacobian_pre) {
868d950fb63SShri Abhyankar       ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr);
869d950fb63SShri Abhyankar     }
870d950fb63SShri Abhyankar 
871d950fb63SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
872d950fb63SShri Abhyankar     snes->linear_its += lits;
873d950fb63SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
874d950fb63SShri Abhyankar     /*
875d950fb63SShri Abhyankar     if (vi->precheckstep) {
876d950fb63SShri Abhyankar       PetscBool changed_y = PETSC_FALSE;
877d950fb63SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
878d950fb63SShri Abhyankar     }
879d950fb63SShri Abhyankar 
880d950fb63SShri Abhyankar     if (PetscLogPrintInfo){
881d950fb63SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
882d950fb63SShri Abhyankar     }
883d950fb63SShri Abhyankar     */
884d950fb63SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
885d950fb63SShri Abhyankar          Y <- X - lambda*Y
886d950fb63SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
887d950fb63SShri Abhyankar     */
888d950fb63SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
889fe835674SShri Abhyankar     ynorm = 1; gnorm = fnorm;
890fe835674SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
891fe835674SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
8922b4ed282SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
8932b4ed282SShri Abhyankar     if (snes->domainerror) {
8942b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
8952b4ed282SShri Abhyankar       PetscFunctionReturn(0);
8962b4ed282SShri Abhyankar     }
8972b4ed282SShri Abhyankar     if (!lssucceed) {
8982b4ed282SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
8992b4ed282SShri Abhyankar 	PetscBool ismin;
9002b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
9012b4ed282SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
9022b4ed282SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
9032b4ed282SShri Abhyankar         break;
9042b4ed282SShri Abhyankar       }
9052b4ed282SShri Abhyankar     }
9062b4ed282SShri Abhyankar     /* Update function and solution vectors */
907fe835674SShri Abhyankar     fnorm = gnorm;
908fe835674SShri Abhyankar     ierr = VecCopy(G,F);CHKERRQ(ierr);
9092b4ed282SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
9102b4ed282SShri Abhyankar     /* Monitor convergence */
9112b4ed282SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
9122b4ed282SShri Abhyankar     snes->iter = i+1;
913fe835674SShri Abhyankar     snes->norm = fnorm;
9142b4ed282SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
9152b4ed282SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
9162b4ed282SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
9172b4ed282SShri Abhyankar     /* Test for convergence, xnorm = || X || */
9182b4ed282SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
919fe835674SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
9202b4ed282SShri Abhyankar     if (snes->reason) break;
9212b4ed282SShri Abhyankar   }
9222b4ed282SShri Abhyankar   if (i == maxits) {
9232b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
9242b4ed282SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
9252b4ed282SShri Abhyankar   }
9262b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9272b4ed282SShri Abhyankar }
9282b4ed282SShri Abhyankar 
9292b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
9302b4ed282SShri Abhyankar /*
9312b4ed282SShri Abhyankar    SNESSetUp_VI - Sets up the internal data structures for the later use
9322b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
9332b4ed282SShri Abhyankar 
9342b4ed282SShri Abhyankar    Input Parameter:
9352b4ed282SShri Abhyankar .  snes - the SNES context
9362b4ed282SShri Abhyankar .  x - the solution vector
9372b4ed282SShri Abhyankar 
9382b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
9392b4ed282SShri Abhyankar 
9402b4ed282SShri Abhyankar    Notes:
9412b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
9422b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
9432b4ed282SShri Abhyankar    the call to SNESSolve().
9442b4ed282SShri Abhyankar  */
9452b4ed282SShri Abhyankar #undef __FUNCT__
9462b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
9472b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
9482b4ed282SShri Abhyankar {
9492b4ed282SShri Abhyankar   PetscErrorCode ierr;
9502b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*) snes->data;
9512b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
9522b4ed282SShri Abhyankar 
9532b4ed282SShri Abhyankar   PetscFunctionBegin;
9542b4ed282SShri Abhyankar   if (!snes->vec_sol_update) {
9552b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
9562b4ed282SShri Abhyankar     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
9572b4ed282SShri Abhyankar   }
9582b4ed282SShri Abhyankar   if (!snes->work) {
9592b4ed282SShri Abhyankar     snes->nwork = 3;
9602b4ed282SShri Abhyankar     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
9612b4ed282SShri Abhyankar     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
9622b4ed282SShri Abhyankar   }
9632b4ed282SShri Abhyankar 
9642b4ed282SShri Abhyankar   /* If the lower and upper bound on variables are not set, set it to
9652b4ed282SShri Abhyankar      -Inf and Inf */
9662b4ed282SShri Abhyankar   if (!vi->xl && !vi->xu) {
9672b4ed282SShri Abhyankar     vi->usersetxbounds = PETSC_FALSE;
9682b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr);
9692b4ed282SShri Abhyankar     ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr);
9702b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr);
9712b4ed282SShri Abhyankar     ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr);
9722b4ed282SShri Abhyankar   } else {
9732b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
9742b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
9752b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
9762b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
9772b4ed282SShri 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]))
9782b4ed282SShri 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.");
9792b4ed282SShri Abhyankar   }
9802b4ed282SShri Abhyankar 
981fe835674SShri Abhyankar   if (snes->ops->solve != SNESSolveVI_RS) {
982fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
983fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr);
984fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr);
985fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr);
986fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
987fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr);
988fe835674SShri Abhyankar 
989fe835674SShri Abhyankar   }
9902b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9912b4ed282SShri Abhyankar }
9922b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
9932b4ed282SShri Abhyankar /*
9942b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
9952b4ed282SShri Abhyankar    with SNESCreate_VI().
9962b4ed282SShri Abhyankar 
9972b4ed282SShri Abhyankar    Input Parameter:
9982b4ed282SShri Abhyankar .  snes - the SNES context
9992b4ed282SShri Abhyankar 
10002b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
10012b4ed282SShri Abhyankar  */
10022b4ed282SShri Abhyankar #undef __FUNCT__
10032b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
10042b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
10052b4ed282SShri Abhyankar {
10062b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*) snes->data;
10072b4ed282SShri Abhyankar   PetscErrorCode ierr;
10082b4ed282SShri Abhyankar 
10092b4ed282SShri Abhyankar   PetscFunctionBegin;
10102b4ed282SShri Abhyankar   if (snes->vec_sol_update) {
10112b4ed282SShri Abhyankar     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
10122b4ed282SShri Abhyankar     snes->vec_sol_update = PETSC_NULL;
10132b4ed282SShri Abhyankar   }
10142b4ed282SShri Abhyankar   if (snes->nwork) {
10152b4ed282SShri Abhyankar     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
10162b4ed282SShri Abhyankar     snes->nwork = 0;
10172b4ed282SShri Abhyankar     snes->work  = PETSC_NULL;
10182b4ed282SShri Abhyankar   }
10192b4ed282SShri Abhyankar 
1020fe835674SShri Abhyankar   if (snes->ops->solve != SNESSolveVI_RS) {
10212b4ed282SShri Abhyankar     /* clear vectors */
1022ee657d29SShri Abhyankar     ierr = VecDestroy(vi->dpsi);CHKERRQ(ierr);
10232b4ed282SShri Abhyankar     ierr = VecDestroy(vi->phi); CHKERRQ(ierr);
10242b4ed282SShri Abhyankar     ierr = VecDestroy(vi->Da); CHKERRQ(ierr);
10252b4ed282SShri Abhyankar     ierr = VecDestroy(vi->Db); CHKERRQ(ierr);
10262b4ed282SShri Abhyankar     ierr = VecDestroy(vi->z); CHKERRQ(ierr);
10272b4ed282SShri Abhyankar     ierr = VecDestroy(vi->t); CHKERRQ(ierr);
10287fe79bc4SShri Abhyankar   }
10297fe79bc4SShri Abhyankar 
10302b4ed282SShri Abhyankar   if (!vi->usersetxbounds) {
10312b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xl); CHKERRQ(ierr);
10322b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xu); CHKERRQ(ierr);
10332b4ed282SShri Abhyankar   }
10347fe79bc4SShri Abhyankar 
1035c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1036c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
10372b4ed282SShri Abhyankar   }
10382b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
10392b4ed282SShri Abhyankar 
10402b4ed282SShri Abhyankar   /* clear composed functions */
10412b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
1042c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
10432b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10442b4ed282SShri Abhyankar }
10457fe79bc4SShri Abhyankar 
10462b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10472b4ed282SShri Abhyankar #undef __FUNCT__
10482b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
10492b4ed282SShri Abhyankar 
10502b4ed282SShri Abhyankar /*
10517fe79bc4SShri Abhyankar   This routine does not actually do a line search but it takes a full newton
10527fe79bc4SShri Abhyankar   step while ensuring that the new iterates remain within the constraints.
10532b4ed282SShri Abhyankar 
10542b4ed282SShri Abhyankar */
1055ace3abfcSBarry 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)
10562b4ed282SShri Abhyankar {
10572b4ed282SShri Abhyankar   PetscErrorCode ierr;
10582b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1059ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
10602b4ed282SShri Abhyankar 
10612b4ed282SShri Abhyankar   PetscFunctionBegin;
10622b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
10632b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
10642b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
10652b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1066c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1067c1a5e521SShri Abhyankar   if (vi->postcheckstep) {
1068c1a5e521SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1069c1a5e521SShri Abhyankar   }
1070c1a5e521SShri Abhyankar   if (changed_y) {
1071c1a5e521SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
10727fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1073c1a5e521SShri Abhyankar   }
1074c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1075c1a5e521SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1076c1a5e521SShri Abhyankar   if (!snes->domainerror) {
10777fe79bc4SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_RS) {
10787fe79bc4SShri Abhyankar        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
10797fe79bc4SShri Abhyankar     } else {
1080c1a5e521SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
10817fe79bc4SShri Abhyankar     }
1082c1a5e521SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1083c1a5e521SShri Abhyankar   }
1084c1a5e521SShri Abhyankar   if (vi->lsmonitor) {
1085c1a5e521SShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1086c1a5e521SShri Abhyankar   }
1087c1a5e521SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1088c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
1089c1a5e521SShri Abhyankar }
1090c1a5e521SShri Abhyankar 
1091c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
1092c1a5e521SShri Abhyankar #undef __FUNCT__
10932b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
10942b4ed282SShri Abhyankar 
10952b4ed282SShri Abhyankar /*
10962b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
10972b4ed282SShri Abhyankar */
1098ace3abfcSBarry 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)
10992b4ed282SShri Abhyankar {
11002b4ed282SShri Abhyankar   PetscErrorCode ierr;
11012b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1102ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
11032b4ed282SShri Abhyankar 
11042b4ed282SShri Abhyankar   PetscFunctionBegin;
11052b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
11062b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11072b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
11087fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
11092b4ed282SShri Abhyankar   if (vi->postcheckstep) {
11102b4ed282SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
11112b4ed282SShri Abhyankar   }
11122b4ed282SShri Abhyankar   if (changed_y) {
11132b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
11147fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
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 }
11247fe79bc4SShri Abhyankar 
11252b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
11262b4ed282SShri Abhyankar #undef __FUNCT__
11272b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
11282b4ed282SShri Abhyankar /*
11297fe79bc4SShri Abhyankar   This routine implements a cubic line search while doing a projection on the variable bounds
11302b4ed282SShri Abhyankar */
1131ace3abfcSBarry 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)
11322b4ed282SShri Abhyankar {
1133fe835674SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
1134fe835674SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
1135fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
1136fe835674SShri Abhyankar   PetscScalar    cinitslope;
1137fe835674SShri Abhyankar #endif
1138fe835674SShri Abhyankar   PetscErrorCode ierr;
1139fe835674SShri Abhyankar   PetscInt       count;
1140fe835674SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1141fe835674SShri Abhyankar   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1142fe835674SShri Abhyankar   MPI_Comm       comm;
1143fe835674SShri Abhyankar 
1144fe835674SShri Abhyankar   PetscFunctionBegin;
1145fe835674SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1146fe835674SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1147fe835674SShri Abhyankar   *flag   = PETSC_TRUE;
1148fe835674SShri Abhyankar 
1149fe835674SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1150fe835674SShri Abhyankar   if (*ynorm == 0.0) {
1151fe835674SShri Abhyankar     if (vi->lsmonitor) {
1152fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
1153fe835674SShri Abhyankar     }
1154fe835674SShri Abhyankar     *gnorm = fnorm;
1155fe835674SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1156fe835674SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1157fe835674SShri Abhyankar     *flag  = PETSC_FALSE;
1158fe835674SShri Abhyankar     goto theend1;
1159fe835674SShri Abhyankar   }
1160fe835674SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1161fe835674SShri Abhyankar     if (vi->lsmonitor) {
1162fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
1163fe835674SShri Abhyankar     }
1164fe835674SShri Abhyankar     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1165fe835674SShri Abhyankar     *ynorm = vi->maxstep;
1166fe835674SShri Abhyankar   }
1167fe835674SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1168fe835674SShri Abhyankar   minlambda = vi->minlambda/rellength;
1169fe835674SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1170fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
1171fe835674SShri Abhyankar   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1172fe835674SShri Abhyankar   initslope = PetscRealPart(cinitslope);
1173fe835674SShri Abhyankar #else
1174fe835674SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1175fe835674SShri Abhyankar #endif
1176fe835674SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
1177fe835674SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
1178fe835674SShri Abhyankar 
1179fe835674SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1180fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1181fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
1182fe835674SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1183fe835674SShri Abhyankar     *flag = PETSC_FALSE;
1184fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1185fe835674SShri Abhyankar     goto theend1;
1186fe835674SShri Abhyankar   }
1187fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
11887fe79bc4SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
1189fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
11907fe79bc4SShri Abhyankar   } else {
11917fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
11927fe79bc4SShri Abhyankar   }
1193fe835674SShri Abhyankar   if (snes->domainerror) {
1194fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1195fe835674SShri Abhyankar     PetscFunctionReturn(0);
1196fe835674SShri Abhyankar   }
1197fe835674SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1198fe835674SShri Abhyankar   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1199fe835674SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1200fe835674SShri Abhyankar     if (vi->lsmonitor) {
1201fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1202fe835674SShri Abhyankar     }
1203fe835674SShri Abhyankar     goto theend1;
1204fe835674SShri Abhyankar   }
1205fe835674SShri Abhyankar 
1206fe835674SShri Abhyankar   /* Fit points with quadratic */
1207fe835674SShri Abhyankar   lambda     = 1.0;
1208fe835674SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
1209fe835674SShri Abhyankar   lambdaprev = lambda;
1210fe835674SShri Abhyankar   gnormprev  = *gnorm;
1211fe835674SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1212fe835674SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
1213fe835674SShri Abhyankar   else                         lambda = lambdatemp;
1214fe835674SShri Abhyankar 
1215fe835674SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1216fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1217fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
1218fe835674SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
1219fe835674SShri Abhyankar     *flag = PETSC_FALSE;
1220fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1221fe835674SShri Abhyankar     goto theend1;
1222fe835674SShri Abhyankar   }
1223fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12247fe79bc4SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
1225fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
12267fe79bc4SShri Abhyankar   } else {
12277fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12287fe79bc4SShri Abhyankar   }
1229fe835674SShri Abhyankar   if (snes->domainerror) {
1230fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1231fe835674SShri Abhyankar     PetscFunctionReturn(0);
1232fe835674SShri Abhyankar   }
1233fe835674SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1234fe835674SShri Abhyankar   if (vi->lsmonitor) {
1235fe835674SShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
1236fe835674SShri Abhyankar   }
1237fe835674SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1238fe835674SShri Abhyankar     if (vi->lsmonitor) {
1239fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
1240fe835674SShri Abhyankar     }
1241fe835674SShri Abhyankar     goto theend1;
1242fe835674SShri Abhyankar   }
1243fe835674SShri Abhyankar 
1244fe835674SShri Abhyankar   /* Fit points with cubic */
1245fe835674SShri Abhyankar   count = 1;
1246fe835674SShri Abhyankar   while (PETSC_TRUE) {
1247fe835674SShri Abhyankar     if (lambda <= minlambda) {
1248fe835674SShri Abhyankar       if (vi->lsmonitor) {
1249fe835674SShri Abhyankar 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
1250fe835674SShri 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);
1251fe835674SShri Abhyankar       }
1252fe835674SShri Abhyankar       *flag = PETSC_FALSE;
1253fe835674SShri Abhyankar       break;
1254fe835674SShri Abhyankar     }
1255fe835674SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
1256fe835674SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
1257fe835674SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1258fe835674SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1259fe835674SShri Abhyankar     d  = b*b - 3*a*initslope;
1260fe835674SShri Abhyankar     if (d < 0.0) d = 0.0;
1261fe835674SShri Abhyankar     if (a == 0.0) {
1262fe835674SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
1263fe835674SShri Abhyankar     } else {
1264fe835674SShri Abhyankar       lambdatemp = (-b + sqrt(d))/(3.0*a);
1265fe835674SShri Abhyankar     }
1266fe835674SShri Abhyankar     lambdaprev = lambda;
1267fe835674SShri Abhyankar     gnormprev  = *gnorm;
1268fe835674SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1269fe835674SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
1270fe835674SShri Abhyankar     else                         lambda     = lambdatemp;
1271fe835674SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1272fe835674SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1273fe835674SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
1274fe835674SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
1275fe835674SShri 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);
1276fe835674SShri Abhyankar       *flag = PETSC_FALSE;
1277fe835674SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1278fe835674SShri Abhyankar       break;
1279fe835674SShri Abhyankar     }
1280fe835674SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12817fe79bc4SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_RS) {
1282fe835674SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
12837fe79bc4SShri Abhyankar     } else {
12847fe79bc4SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12857fe79bc4SShri Abhyankar     }
1286fe835674SShri Abhyankar     if (snes->domainerror) {
1287fe835674SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1288fe835674SShri Abhyankar       PetscFunctionReturn(0);
1289fe835674SShri Abhyankar     }
1290fe835674SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1291fe835674SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
1292fe835674SShri Abhyankar       if (vi->lsmonitor) {
1293fe835674SShri Abhyankar 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1294fe835674SShri Abhyankar       }
1295fe835674SShri Abhyankar       break;
1296fe835674SShri Abhyankar     } else {
1297fe835674SShri Abhyankar       if (vi->lsmonitor) {
1298fe835674SShri Abhyankar         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1299fe835674SShri Abhyankar       }
1300fe835674SShri Abhyankar     }
1301fe835674SShri Abhyankar     count++;
1302fe835674SShri Abhyankar   }
1303fe835674SShri Abhyankar   theend1:
1304fe835674SShri Abhyankar   /* Optional user-defined check for line search step validity */
1305fe835674SShri Abhyankar   if (vi->postcheckstep && *flag) {
1306fe835674SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1307fe835674SShri Abhyankar     if (changed_y) {
1308fe835674SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1309fe835674SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1310fe835674SShri Abhyankar     }
1311fe835674SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
1312fe835674SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13137fe79bc4SShri Abhyankar       if (snes->ops->solve == SNESSolveVI_RS) {
1314fe835674SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
13157fe79bc4SShri Abhyankar       } else {
13167fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13177fe79bc4SShri Abhyankar       }
1318fe835674SShri Abhyankar       if (snes->domainerror) {
1319fe835674SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1320fe835674SShri Abhyankar         PetscFunctionReturn(0);
1321fe835674SShri Abhyankar       }
1322fe835674SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1323fe835674SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
1324fe835674SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
1325fe835674SShri Abhyankar     }
1326fe835674SShri Abhyankar   }
1327fe835674SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1328fe835674SShri Abhyankar   PetscFunctionReturn(0);
1329fe835674SShri Abhyankar }
1330fe835674SShri Abhyankar 
13312b4ed282SShri Abhyankar #undef __FUNCT__
13322b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
13332b4ed282SShri Abhyankar /*
13347fe79bc4SShri Abhyankar   This routine does a quadratic line search while keeping the iterates within the variable bounds
13352b4ed282SShri Abhyankar */
1336ace3abfcSBarry 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)
13372b4ed282SShri Abhyankar {
13382b4ed282SShri Abhyankar   /*
13392b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
13402b4ed282SShri Abhyankar      minimization problem:
13412b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
13422b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
13432b4ed282SShri Abhyankar    */
13442b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
13452b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
13462b4ed282SShri Abhyankar   PetscScalar    cinitslope;
13472b4ed282SShri Abhyankar #endif
13482b4ed282SShri Abhyankar   PetscErrorCode ierr;
13492b4ed282SShri Abhyankar   PetscInt       count;
13502b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1351ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
13522b4ed282SShri Abhyankar 
13532b4ed282SShri Abhyankar   PetscFunctionBegin;
13542b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
13552b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
13562b4ed282SShri Abhyankar 
13572b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
13582b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
1359c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1360c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
13612b4ed282SShri Abhyankar     }
13622b4ed282SShri Abhyankar     *gnorm = fnorm;
13632b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
13642b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
13652b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
13662b4ed282SShri Abhyankar     goto theend2;
13672b4ed282SShri Abhyankar   }
13682b4ed282SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
13692b4ed282SShri Abhyankar     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
13702b4ed282SShri Abhyankar     *ynorm = vi->maxstep;
13712b4ed282SShri Abhyankar   }
13722b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
13732b4ed282SShri Abhyankar   minlambda = vi->minlambda/rellength;
13747fe79bc4SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
13752b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
13767fe79bc4SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
13772b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
13782b4ed282SShri Abhyankar #else
13797fe79bc4SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
13802b4ed282SShri Abhyankar #endif
13812b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
13822b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
13832b4ed282SShri Abhyankar   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
13842b4ed282SShri Abhyankar 
13852b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
13867fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
13872b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
13882b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
13892b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
13902b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
13912b4ed282SShri Abhyankar     goto theend2;
13922b4ed282SShri Abhyankar   }
13932b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13947fe79bc4SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
13957fe79bc4SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
13967fe79bc4SShri Abhyankar   } else {
13977fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13987fe79bc4SShri Abhyankar   }
13992b4ed282SShri Abhyankar   if (snes->domainerror) {
14002b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14012b4ed282SShri Abhyankar     PetscFunctionReturn(0);
14022b4ed282SShri Abhyankar   }
14032b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
14042b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1405c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1406c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
14072b4ed282SShri Abhyankar     }
14082b4ed282SShri Abhyankar     goto theend2;
14092b4ed282SShri Abhyankar   }
14102b4ed282SShri Abhyankar 
14112b4ed282SShri Abhyankar   /* Fit points with quadratic */
14122b4ed282SShri Abhyankar   lambda = 1.0;
14132b4ed282SShri Abhyankar   count = 1;
14142b4ed282SShri Abhyankar   while (PETSC_TRUE) {
14152b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
1416c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1417c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
1418c92abb8aSShri 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);
14192b4ed282SShri Abhyankar       }
14202b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
14212b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
14222b4ed282SShri Abhyankar       break;
14232b4ed282SShri Abhyankar     }
14242b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
14252b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
14262b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
14272b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
14282b4ed282SShri Abhyankar 
14292b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
14307fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
14312b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
14322b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
14332b4ed282SShri 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);
14342b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
14352b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
14362b4ed282SShri Abhyankar       break;
14372b4ed282SShri Abhyankar     }
14382b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
14392b4ed282SShri Abhyankar     if (snes->domainerror) {
14402b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14412b4ed282SShri Abhyankar       PetscFunctionReturn(0);
14422b4ed282SShri Abhyankar     }
14437fe79bc4SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_RS) {
14447fe79bc4SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
14457fe79bc4SShri Abhyankar     } else {
14462b4ed282SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
14477fe79bc4SShri Abhyankar     }
14482b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
14492b4ed282SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1450c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1451c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
14522b4ed282SShri Abhyankar       }
14532b4ed282SShri Abhyankar       break;
14542b4ed282SShri Abhyankar     }
14552b4ed282SShri Abhyankar     count++;
14562b4ed282SShri Abhyankar   }
14572b4ed282SShri Abhyankar   theend2:
14582b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
14592b4ed282SShri Abhyankar   if (vi->postcheckstep) {
14602b4ed282SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
14612b4ed282SShri Abhyankar     if (changed_y) {
14622b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
14637fe79bc4SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
14642b4ed282SShri Abhyankar     }
14652b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
14662b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
14672b4ed282SShri Abhyankar       if (snes->domainerror) {
14682b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14692b4ed282SShri Abhyankar         PetscFunctionReturn(0);
14702b4ed282SShri Abhyankar       }
14717fe79bc4SShri Abhyankar       if (snes->ops->solve == SNESSolveVI_RS) {
14727fe79bc4SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
14737fe79bc4SShri Abhyankar       } else {
14747fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
14757fe79bc4SShri Abhyankar       }
14767fe79bc4SShri Abhyankar 
14772b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
14782b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
14792b4ed282SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
14802b4ed282SShri Abhyankar     }
14812b4ed282SShri Abhyankar   }
14822b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14832b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14842b4ed282SShri Abhyankar }
14852b4ed282SShri Abhyankar 
1486ace3abfcSBarry 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*/
14872b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
14882b4ed282SShri Abhyankar EXTERN_C_BEGIN
14892b4ed282SShri Abhyankar #undef __FUNCT__
14902b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI"
14917087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
14922b4ed282SShri Abhyankar {
14932b4ed282SShri Abhyankar   PetscFunctionBegin;
14942b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->LineSearch = func;
14952b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->lsP        = lsctx;
14962b4ed282SShri Abhyankar   PetscFunctionReturn(0);
14972b4ed282SShri Abhyankar }
14982b4ed282SShri Abhyankar EXTERN_C_END
14992b4ed282SShri Abhyankar 
15002b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15012b4ed282SShri Abhyankar EXTERN_C_BEGIN
15022b4ed282SShri Abhyankar #undef __FUNCT__
15032b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
15047087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
15052b4ed282SShri Abhyankar {
15062b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
15072b4ed282SShri Abhyankar   PetscErrorCode ierr;
15082b4ed282SShri Abhyankar 
15092b4ed282SShri Abhyankar   PetscFunctionBegin;
1510c92abb8aSShri Abhyankar   if (flg && !vi->lsmonitor) {
1511c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr);
1512c92abb8aSShri Abhyankar   } else if (!flg && vi->lsmonitor) {
1513c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
15142b4ed282SShri Abhyankar   }
15152b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15162b4ed282SShri Abhyankar }
15172b4ed282SShri Abhyankar EXTERN_C_END
15182b4ed282SShri Abhyankar 
15192b4ed282SShri Abhyankar /*
15202b4ed282SShri Abhyankar    SNESView_VI - Prints info from the SNESVI data structure.
15212b4ed282SShri Abhyankar 
15222b4ed282SShri Abhyankar    Input Parameters:
15232b4ed282SShri Abhyankar .  SNES - the SNES context
15242b4ed282SShri Abhyankar .  viewer - visualization context
15252b4ed282SShri Abhyankar 
15262b4ed282SShri Abhyankar    Application Interface Routine: SNESView()
15272b4ed282SShri Abhyankar */
15282b4ed282SShri Abhyankar #undef __FUNCT__
15292b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI"
15302b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
15312b4ed282SShri Abhyankar {
15322b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
153378c4b9d3SShri Abhyankar   const char     *cstr,*tstr;
15342b4ed282SShri Abhyankar   PetscErrorCode ierr;
1535ace3abfcSBarry Smith   PetscBool     iascii;
15362b4ed282SShri Abhyankar 
15372b4ed282SShri Abhyankar   PetscFunctionBegin;
15382b4ed282SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
15392b4ed282SShri Abhyankar   if (iascii) {
15402b4ed282SShri Abhyankar     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
15412b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
15422b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
15432b4ed282SShri Abhyankar     else                                                             cstr = "unknown";
154478c4b9d3SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_SS)      tstr = "Semismooth";
15450a7c7af0SShri Abhyankar     else if (snes->ops->solve == SNESSolveVI_RS) tstr = "Reduced Space";
154678c4b9d3SShri Abhyankar     else                                         tstr = "unknown";
154778c4b9d3SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  VI algorithm: %s\n",tstr);CHKERRQ(ierr);
15482b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
15492b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
15502b4ed282SShri Abhyankar   } else {
15512b4ed282SShri Abhyankar     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
15522b4ed282SShri Abhyankar   }
15532b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15542b4ed282SShri Abhyankar }
15552b4ed282SShri Abhyankar 
15565559b345SBarry Smith #undef __FUNCT__
15575559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
15585559b345SBarry Smith /*@
15592b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
15602b4ed282SShri Abhyankar 
15612b4ed282SShri Abhyankar    Input Parameters:
15622b4ed282SShri Abhyankar .  snes - the SNES context.
15632b4ed282SShri Abhyankar .  xl   - lower bound.
15642b4ed282SShri Abhyankar .  xu   - upper bound.
15652b4ed282SShri Abhyankar 
15662b4ed282SShri Abhyankar    Notes:
15672b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
156884c105d7SBarry Smith    PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp().
156984c105d7SBarry Smith 
15705559b345SBarry Smith @*/
157169c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
15722b4ed282SShri Abhyankar {
15732b4ed282SShri Abhyankar   SNES_VI  *vi = (SNES_VI*)snes->data;
15742b4ed282SShri Abhyankar 
15752b4ed282SShri Abhyankar   PetscFunctionBegin;
15762b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
15772b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
15782b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
15792b4ed282SShri Abhyankar   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
15802b4ed282SShri 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);
15812b4ed282SShri 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);
15825559b345SBarry Smith 
15832b4ed282SShri Abhyankar   vi->usersetxbounds = PETSC_TRUE;
15842b4ed282SShri Abhyankar   vi->xl = xl;
15852b4ed282SShri Abhyankar   vi->xu = xu;
15862b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15872b4ed282SShri Abhyankar }
15882b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15892b4ed282SShri Abhyankar /*
15902b4ed282SShri Abhyankar    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
15912b4ed282SShri Abhyankar 
15922b4ed282SShri Abhyankar    Input Parameter:
15932b4ed282SShri Abhyankar .  snes - the SNES context
15942b4ed282SShri Abhyankar 
15952b4ed282SShri Abhyankar    Application Interface Routine: SNESSetFromOptions()
15962b4ed282SShri Abhyankar */
15972b4ed282SShri Abhyankar #undef __FUNCT__
15982b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI"
15992b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
16002b4ed282SShri Abhyankar {
16012b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
16027fe79bc4SShri Abhyankar   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
16035559b345SBarry Smith   const char     *vies[] = {"ss","rs"};
16042b4ed282SShri Abhyankar   PetscErrorCode ierr;
16052b4ed282SShri Abhyankar   PetscInt       indx;
160669c03620SShri Abhyankar   PetscBool      flg,set,flg2;
16072b4ed282SShri Abhyankar 
16082b4ed282SShri Abhyankar   PetscFunctionBegin;
16092b4ed282SShri Abhyankar   ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
16109308c137SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
16119308c137SBarry Smith   if (flg) {
16129308c137SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
16139308c137SBarry Smith   }
16142b4ed282SShri Abhyankar   ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
16152b4ed282SShri Abhyankar   ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
16162b4ed282SShri Abhyankar   ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
16172b4ed282SShri Abhyankar   ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
1618acfcf0e5SJed Brown   ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
16192b4ed282SShri Abhyankar   if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
16205559b345SBarry Smith   ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr);
162169c03620SShri Abhyankar   if (flg2) {
162269c03620SShri Abhyankar     switch (indx) {
162369c03620SShri Abhyankar     case 0:
162469c03620SShri Abhyankar       snes->ops->solve = SNESSolveVI_SS;
162569c03620SShri Abhyankar       break;
162669c03620SShri Abhyankar     case 1:
1627d950fb63SShri Abhyankar       snes->ops->solve = SNESSolveVI_RS;
1628d950fb63SShri Abhyankar       break;
162969c03620SShri Abhyankar     }
163069c03620SShri Abhyankar   }
16317fe79bc4SShri Abhyankar   ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"basic",&indx,&flg);CHKERRQ(ierr);
16322b4ed282SShri Abhyankar   if (flg) {
16332b4ed282SShri Abhyankar     switch (indx) {
16342b4ed282SShri Abhyankar     case 0:
16352b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
16362b4ed282SShri Abhyankar       break;
16372b4ed282SShri Abhyankar     case 1:
16382b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
16392b4ed282SShri Abhyankar       break;
16402b4ed282SShri Abhyankar     case 2:
16412b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
16422b4ed282SShri Abhyankar       break;
16432b4ed282SShri Abhyankar     case 3:
16442b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
16452b4ed282SShri Abhyankar       break;
16462b4ed282SShri Abhyankar     }
1647fe835674SShri Abhyankar   }
16482b4ed282SShri Abhyankar   ierr = PetscOptionsTail();CHKERRQ(ierr);
16492b4ed282SShri Abhyankar   PetscFunctionReturn(0);
16502b4ed282SShri Abhyankar }
16512b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
16522b4ed282SShri Abhyankar /*MC
16532b4ed282SShri Abhyankar       SNESVI - Semismooth newton method based nonlinear solver that uses a line search
16542b4ed282SShri Abhyankar 
16552b4ed282SShri Abhyankar    Options Database:
16562b4ed282SShri Abhyankar +   -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search
16572b4ed282SShri Abhyankar .   -snes_vi_alpha <alpha> - Sets alpha
16582b4ed282SShri 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)
16592b4ed282SShri Abhyankar .   -snes_vi_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
16602b4ed282SShri Abhyankar -   -snes_vi_monitor - print information about progress of line searches
16612b4ed282SShri Abhyankar 
16622b4ed282SShri Abhyankar 
16632b4ed282SShri Abhyankar    Level: beginner
16642b4ed282SShri Abhyankar 
16652b4ed282SShri Abhyankar .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
16662b4ed282SShri Abhyankar            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
16672b4ed282SShri Abhyankar            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
16682b4ed282SShri Abhyankar 
16692b4ed282SShri Abhyankar M*/
16702b4ed282SShri Abhyankar EXTERN_C_BEGIN
16712b4ed282SShri Abhyankar #undef __FUNCT__
16722b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI"
16737087cfbeSBarry Smith PetscErrorCode  SNESCreate_VI(SNES snes)
16742b4ed282SShri Abhyankar {
16752b4ed282SShri Abhyankar   PetscErrorCode ierr;
16762b4ed282SShri Abhyankar   SNES_VI      *vi;
16772b4ed282SShri Abhyankar 
16782b4ed282SShri Abhyankar   PetscFunctionBegin;
16792b4ed282SShri Abhyankar   snes->ops->setup	     = SNESSetUp_VI;
168069c03620SShri Abhyankar   snes->ops->solve	     = SNESSolveVI_SS;
16812b4ed282SShri Abhyankar   snes->ops->destroy	     = SNESDestroy_VI;
16822b4ed282SShri Abhyankar   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
16832b4ed282SShri Abhyankar   snes->ops->view            = SNESView_VI;
16842b4ed282SShri Abhyankar   snes->ops->converged       = SNESDefaultConverged_VI;
16852b4ed282SShri Abhyankar 
16862b4ed282SShri Abhyankar   ierr                   = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
16872b4ed282SShri Abhyankar   snes->data    	 = (void*)vi;
16882b4ed282SShri Abhyankar   vi->alpha		 = 1.e-4;
16892b4ed282SShri Abhyankar   vi->maxstep		 = 1.e8;
16902b4ed282SShri Abhyankar   vi->minlambda         = 1.e-12;
16917fe79bc4SShri Abhyankar   vi->LineSearch        = SNESLineSearchNo_VI;
16922b4ed282SShri Abhyankar   vi->lsP               = PETSC_NULL;
16932b4ed282SShri Abhyankar   vi->postcheckstep     = PETSC_NULL;
16942b4ed282SShri Abhyankar   vi->postcheck         = PETSC_NULL;
16952b4ed282SShri Abhyankar   vi->precheckstep      = PETSC_NULL;
16962b4ed282SShri Abhyankar   vi->precheck          = PETSC_NULL;
16972b4ed282SShri Abhyankar   vi->const_tol         =  2.2204460492503131e-16;
16982b4ed282SShri Abhyankar 
16992b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
17002b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
17012b4ed282SShri Abhyankar 
17022b4ed282SShri Abhyankar   PetscFunctionReturn(0);
17032b4ed282SShri Abhyankar }
17042b4ed282SShri Abhyankar EXTERN_C_END
1705