xref: /petsc/src/snes/impls/vi/vi.c (revision 1130f37da7db8287cb85ef5e2737a21840517f1b)
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"
527d4218bSShri Abhyankar #include "../include/private/matimpl.h"
62b4ed282SShri Abhyankar 
79308c137SBarry Smith #undef __FUNCT__
89308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
97087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
109308c137SBarry Smith {
119308c137SBarry Smith   PetscErrorCode          ierr;
129308c137SBarry Smith   SNES_VI                 *vi = (SNES_VI*)snes->data;
139308c137SBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;
149308c137SBarry Smith   const PetscScalar       *x,*xl,*xu,*f;
1509573ac7SBarry Smith   PetscInt                i,n,act = 0;
169308c137SBarry Smith   PetscReal               rnorm,fnorm;
179308c137SBarry Smith 
189308c137SBarry Smith   PetscFunctionBegin;
199308c137SBarry Smith   if (!dummy) {
209308c137SBarry Smith     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
219308c137SBarry Smith   }
229308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
239308c137SBarry Smith   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
249308c137SBarry Smith   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
259308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
263f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
279308c137SBarry Smith 
289308c137SBarry Smith   rnorm = 0.0;
299308c137SBarry Smith   for (i=0; i<n; i++) {
307b2c10efSShri 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];
3109573ac7SBarry Smith     else act++;
329308c137SBarry Smith   }
333f731af5SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
349308c137SBarry Smith   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
359308c137SBarry Smith   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
369308c137SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
379308c137SBarry Smith   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
389308c137SBarry Smith   fnorm = sqrt(fnorm);
3909573ac7SBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES VI Function norm %14.12e Active constraints %D\n",its,fnorm,act);CHKERRQ(ierr);
409308c137SBarry Smith   if (!dummy) {
419308c137SBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
429308c137SBarry Smith   }
439308c137SBarry Smith   PetscFunctionReturn(0);
449308c137SBarry Smith }
459308c137SBarry Smith 
462b4ed282SShri Abhyankar /*
472b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
482b4ed282SShri 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
492b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
502b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
512b4ed282SShri Abhyankar */
522b4ed282SShri Abhyankar #undef __FUNCT__
532b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
54ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
552b4ed282SShri Abhyankar {
562b4ed282SShri Abhyankar   PetscReal      a1;
572b4ed282SShri Abhyankar   PetscErrorCode ierr;
58ace3abfcSBarry Smith   PetscBool     hastranspose;
592b4ed282SShri Abhyankar 
602b4ed282SShri Abhyankar   PetscFunctionBegin;
612b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
622b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
632b4ed282SShri Abhyankar   if (hastranspose) {
642b4ed282SShri Abhyankar     /* Compute || J^T F|| */
652b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
662b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
672b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
682b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
692b4ed282SShri Abhyankar   } else {
702b4ed282SShri Abhyankar     Vec         work;
712b4ed282SShri Abhyankar     PetscScalar result;
722b4ed282SShri Abhyankar     PetscReal   wnorm;
732b4ed282SShri Abhyankar 
742b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
752b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
762b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
772b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
782b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
792b4ed282SShri Abhyankar     ierr = VecDestroy(work);CHKERRQ(ierr);
802b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
812b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
822b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
832b4ed282SShri Abhyankar   }
842b4ed282SShri Abhyankar   PetscFunctionReturn(0);
852b4ed282SShri Abhyankar }
862b4ed282SShri Abhyankar 
872b4ed282SShri Abhyankar /*
882b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
892b4ed282SShri Abhyankar */
902b4ed282SShri Abhyankar #undef __FUNCT__
912b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
922b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
932b4ed282SShri Abhyankar {
942b4ed282SShri Abhyankar   PetscReal      a1,a2;
952b4ed282SShri Abhyankar   PetscErrorCode ierr;
96ace3abfcSBarry Smith   PetscBool     hastranspose;
972b4ed282SShri Abhyankar 
982b4ed282SShri Abhyankar   PetscFunctionBegin;
992b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1002b4ed282SShri Abhyankar   if (hastranspose) {
1012b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1022b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1032b4ed282SShri Abhyankar 
1042b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1052b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
1062b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
1072b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
1082b4ed282SShri Abhyankar     if (a1 != 0.0) {
1092b4ed282SShri Abhyankar       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
1102b4ed282SShri Abhyankar     }
1112b4ed282SShri Abhyankar   }
1122b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1132b4ed282SShri Abhyankar }
1142b4ed282SShri Abhyankar 
1152b4ed282SShri Abhyankar /*
1162b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
1172b4ed282SShri Abhyankar 
1182b4ed282SShri Abhyankar   Notes:
1192b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1202b4ed282SShri Abhyankar   merit < abstol
1212b4ed282SShri Abhyankar   merit < rtol*merit_initial
1222b4ed282SShri Abhyankar */
1232b4ed282SShri Abhyankar #undef __FUNCT__
1242b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
1257fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
1262b4ed282SShri Abhyankar {
1272b4ed282SShri Abhyankar   PetscErrorCode ierr;
1282b4ed282SShri Abhyankar 
1292b4ed282SShri Abhyankar   PetscFunctionBegin;
1302b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1312b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
1322b4ed282SShri Abhyankar 
1332b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1342b4ed282SShri Abhyankar 
1352b4ed282SShri Abhyankar   if (!it) {
1362b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1377fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
1382b4ed282SShri Abhyankar   }
1397fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
1402b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
1412b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1427fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
1437fe79bc4SShri Abhyankar     ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr);
1442b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
1452b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
1462b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
1472b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
1482b4ed282SShri Abhyankar   }
1492b4ed282SShri Abhyankar 
1502b4ed282SShri Abhyankar   if (it && !*reason) {
1517fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
1527fe79bc4SShri Abhyankar       ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr);
1532b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
1542b4ed282SShri Abhyankar     }
1552b4ed282SShri Abhyankar   }
1562b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1572b4ed282SShri Abhyankar }
1582b4ed282SShri Abhyankar 
1592b4ed282SShri Abhyankar /*
1602b4ed282SShri Abhyankar   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
1612b4ed282SShri Abhyankar 
1622b4ed282SShri Abhyankar   Input Parameter:
1632b4ed282SShri Abhyankar . phi - the semismooth function
1642b4ed282SShri Abhyankar 
1652b4ed282SShri Abhyankar   Output Parameter:
1662b4ed282SShri Abhyankar . merit - the merit function
1672b4ed282SShri Abhyankar . phinorm - ||phi||
1682b4ed282SShri Abhyankar 
1692b4ed282SShri Abhyankar   Notes:
1702b4ed282SShri Abhyankar   The merit function for the mixed complementarity problem is defined as
1712b4ed282SShri Abhyankar      merit = 0.5*phi^T*phi
1722b4ed282SShri Abhyankar */
1732b4ed282SShri Abhyankar #undef __FUNCT__
1742b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction"
1752b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
1762b4ed282SShri Abhyankar {
1772b4ed282SShri Abhyankar   PetscErrorCode ierr;
1782b4ed282SShri Abhyankar 
1792b4ed282SShri Abhyankar   PetscFunctionBegin;
1802b4ed282SShri Abhyankar   ierr = VecNormBegin(phi,NORM_2,phinorm);
1812b4ed282SShri Abhyankar   ierr = VecNormEnd(phi,NORM_2,phinorm);
1822b4ed282SShri Abhyankar 
1832b4ed282SShri Abhyankar   *merit = 0.5*(*phinorm)*(*phinorm);
1842b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1852b4ed282SShri Abhyankar }
1862b4ed282SShri Abhyankar 
1877f6f51e9SSatish Balay PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b)
1884c21c6cdSBarry Smith {
1894c21c6cdSBarry Smith   return a + b - sqrt(a*a + b*b);
1904c21c6cdSBarry Smith }
1914c21c6cdSBarry Smith 
1927f6f51e9SSatish Balay PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b)
1934c21c6cdSBarry Smith {
1945559b345SBarry Smith   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return  1.0 - a/ sqrt(a*a + b*b);
1955559b345SBarry Smith   else return .5;
1964c21c6cdSBarry Smith }
1974c21c6cdSBarry Smith 
1982b4ed282SShri Abhyankar /*
1991f79c099SShri Abhyankar    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
2002b4ed282SShri Abhyankar 
2012b4ed282SShri Abhyankar    Input Parameters:
2022b4ed282SShri Abhyankar .  snes - the SNES context
2032b4ed282SShri Abhyankar .  x - current iterate
2042b4ed282SShri Abhyankar .  functx - user defined function context
2052b4ed282SShri Abhyankar 
2062b4ed282SShri Abhyankar    Output Parameters:
2072b4ed282SShri Abhyankar .  phi - Semismooth function
2082b4ed282SShri Abhyankar 
2092b4ed282SShri Abhyankar */
2102b4ed282SShri Abhyankar #undef __FUNCT__
2111f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction"
2121f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
2132b4ed282SShri Abhyankar {
2142b4ed282SShri Abhyankar   PetscErrorCode  ierr;
2152b4ed282SShri Abhyankar   SNES_VI         *vi = (SNES_VI*)snes->data;
2162b4ed282SShri Abhyankar   Vec             Xl = vi->xl,Xu = vi->xu,F = snes->vec_func;
2174c21c6cdSBarry Smith   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u;
2182b4ed282SShri Abhyankar   PetscInt        i,nlocal;
2192b4ed282SShri Abhyankar 
2202b4ed282SShri Abhyankar   PetscFunctionBegin;
2212b4ed282SShri Abhyankar   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
2222b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2232b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
2242b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
2252b4ed282SShri Abhyankar   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
2262b4ed282SShri Abhyankar   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
2272b4ed282SShri Abhyankar   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
2282b4ed282SShri Abhyankar 
2292b4ed282SShri Abhyankar   for (i=0;i < nlocal;i++) {
2305559b345SBarry Smith     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { /* no constraints on variable */
2314c21c6cdSBarry Smith       phi_arr[i] = f_arr[i];
2325559b345SBarry Smith     } else if (l[i] <= PETSC_VI_NINF) {                      /* upper bound on variable only */
2334c21c6cdSBarry Smith       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
2345559b345SBarry Smith     } else if (u[i] >= PETSC_VI_INF) {                       /* lower bound on variable only */
2354c21c6cdSBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
2365559b345SBarry Smith     } else if (l[i] == u[i]) {
2372b4ed282SShri Abhyankar       phi_arr[i] = l[i] - x_arr[i];
2385559b345SBarry Smith     } else {                                                /* both bounds on variable */
2394c21c6cdSBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
2402b4ed282SShri Abhyankar     }
2412b4ed282SShri Abhyankar   }
2422b4ed282SShri Abhyankar 
2432b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
2442b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
2452b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
2462b4ed282SShri Abhyankar   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
2472b4ed282SShri Abhyankar   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
2482b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2492b4ed282SShri Abhyankar }
2502b4ed282SShri Abhyankar 
2512b4ed282SShri Abhyankar /*
252a79edbefSShri Abhyankar    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
253a79edbefSShri Abhyankar                                           the semismooth jacobian.
2542b4ed282SShri Abhyankar */
2552b4ed282SShri Abhyankar #undef __FUNCT__
256a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
257a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
2582b4ed282SShri Abhyankar {
2592b4ed282SShri Abhyankar   PetscErrorCode ierr;
2602b4ed282SShri Abhyankar   SNES_VI      *vi = (SNES_VI*)snes->data;
2614c21c6cdSBarry Smith   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
2622b4ed282SShri Abhyankar   PetscInt       i,nlocal;
2632b4ed282SShri Abhyankar 
2642b4ed282SShri Abhyankar   PetscFunctionBegin;
2652b4ed282SShri Abhyankar 
2662b4ed282SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
2672b4ed282SShri Abhyankar   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
2682b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr);
2692b4ed282SShri Abhyankar   ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr);
270a79edbefSShri Abhyankar   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
271a79edbefSShri Abhyankar   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
2722b4ed282SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
2734c21c6cdSBarry Smith 
2742b4ed282SShri Abhyankar   for (i=0;i< nlocal;i++) {
2755559b345SBarry Smith     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {/* no constraints on variable */
2764c21c6cdSBarry Smith       da[i] = 0;
2772b4ed282SShri Abhyankar       db[i] = 1;
2785559b345SBarry Smith     } else if (l[i] <= PETSC_VI_NINF) {                     /* upper bound on variable only */
2794c21c6cdSBarry Smith       da[i] = DPhi(u[i] - x[i], -f[i]);
2804c21c6cdSBarry Smith       db[i] = DPhi(-f[i],u[i] - x[i]);
2815559b345SBarry Smith     } else if (u[i] >= PETSC_VI_INF) {                      /* lower bound on variable only */
2825559b345SBarry Smith       da[i] = DPhi(x[i] - l[i], f[i]);
2835559b345SBarry Smith       db[i] = DPhi(f[i],x[i] - l[i]);
2845559b345SBarry Smith     } else if (l[i] == u[i]) {                              /* fixed variable */
2854c21c6cdSBarry Smith       da[i] = 1;
2862b4ed282SShri Abhyankar       db[i] = 0;
2875559b345SBarry Smith     } else {                                                /* upper and lower bounds on variable */
2884c21c6cdSBarry Smith       da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
2894c21c6cdSBarry Smith       db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
2904c21c6cdSBarry Smith       da2 = DPhi(u[i] - x[i], -f[i]);
2914c21c6cdSBarry Smith       db2 = DPhi(-f[i],u[i] - x[i]);
2924c21c6cdSBarry Smith       da[i] = da1 + db1*da2;
2934c21c6cdSBarry Smith       db[i] = db1*db2;
2942b4ed282SShri Abhyankar     }
2952b4ed282SShri Abhyankar   }
2965559b345SBarry Smith 
2972b4ed282SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
2982b4ed282SShri Abhyankar   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
2992b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr);
3002b4ed282SShri Abhyankar   ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr);
301a79edbefSShri Abhyankar   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
302a79edbefSShri Abhyankar   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
303a79edbefSShri Abhyankar   PetscFunctionReturn(0);
304a79edbefSShri Abhyankar }
305a79edbefSShri Abhyankar 
306a79edbefSShri Abhyankar /*
307a79edbefSShri 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.
308a79edbefSShri Abhyankar 
309a79edbefSShri Abhyankar    Input Parameters:
310a79edbefSShri Abhyankar .  Da       - Diagonal shift vector for the semismooth jacobian.
311a79edbefSShri Abhyankar .  Db       - Row scaling vector for the semismooth jacobian.
312a79edbefSShri Abhyankar 
313a79edbefSShri Abhyankar    Output Parameters:
314a79edbefSShri Abhyankar .  jac      - semismooth jacobian
315a79edbefSShri Abhyankar .  jac_pre  - optional preconditioning matrix
316a79edbefSShri Abhyankar 
317a79edbefSShri Abhyankar    Notes:
318a79edbefSShri Abhyankar    The semismooth jacobian matrix is given by
319a79edbefSShri Abhyankar    jac = Da + Db*jacfun
320a79edbefSShri Abhyankar    where Db is the row scaling matrix stored as a vector,
321a79edbefSShri Abhyankar          Da is the diagonal perturbation matrix stored as a vector
322a79edbefSShri Abhyankar    and   jacfun is the jacobian of the original nonlinear function.
323a79edbefSShri Abhyankar */
324a79edbefSShri Abhyankar #undef __FUNCT__
325a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian"
326a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
327a79edbefSShri Abhyankar {
328a79edbefSShri Abhyankar   PetscErrorCode ierr;
329a79edbefSShri Abhyankar 
3302b4ed282SShri Abhyankar   /* Do row scaling  and add diagonal perturbation */
331a79edbefSShri Abhyankar   ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr);
332a79edbefSShri Abhyankar   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
333a79edbefSShri Abhyankar   if (jac != jac_pre) { /* If jac and jac_pre are different */
334a79edbefSShri Abhyankar     ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL);
335a79edbefSShri Abhyankar     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
3362b4ed282SShri Abhyankar   }
3372b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3382b4ed282SShri Abhyankar }
3392b4ed282SShri Abhyankar 
3402b4ed282SShri Abhyankar /*
341ee657d29SShri Abhyankar    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
342ee657d29SShri Abhyankar 
343ee657d29SShri Abhyankar    Input Parameters:
344ee657d29SShri Abhyankar    phi - semismooth function.
345ee657d29SShri Abhyankar    H   - semismooth jacobian
346ee657d29SShri Abhyankar 
347ee657d29SShri Abhyankar    Output Parameters:
348ee657d29SShri Abhyankar    dpsi - merit function gradient
349ee657d29SShri Abhyankar 
350ee657d29SShri Abhyankar    Notes:
351ee657d29SShri Abhyankar   The merit function gradient is computed as follows
352ee657d29SShri Abhyankar         dpsi = H^T*phi
353ee657d29SShri Abhyankar */
354ee657d29SShri Abhyankar #undef __FUNCT__
355ee657d29SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
356ee657d29SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
357ee657d29SShri Abhyankar {
358ee657d29SShri Abhyankar   PetscErrorCode ierr;
359ee657d29SShri Abhyankar 
360ee657d29SShri Abhyankar   PetscFunctionBegin;
361ee657d29SShri Abhyankar   ierr = MatMultTranspose(H,phi,dpsi);
362ee657d29SShri Abhyankar   PetscFunctionReturn(0);
363ee657d29SShri Abhyankar }
364ee657d29SShri Abhyankar 
365c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
366c1a5e521SShri Abhyankar /*
367c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
368c1a5e521SShri Abhyankar 
369c1a5e521SShri Abhyankar    Input Parameters:
370c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
371c1a5e521SShri Abhyankar 
372c1a5e521SShri Abhyankar    Output Parameters:
373c1a5e521SShri Abhyankar .  X - Bound projected X
374c1a5e521SShri Abhyankar 
375c1a5e521SShri Abhyankar */
376c1a5e521SShri Abhyankar 
377c1a5e521SShri Abhyankar #undef __FUNCT__
378c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
379c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
380c1a5e521SShri Abhyankar {
381c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
382c1a5e521SShri Abhyankar   SNES_VI           *vi = (SNES_VI*)snes->data;
383c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
384c1a5e521SShri Abhyankar   PetscScalar       *x;
385c1a5e521SShri Abhyankar   PetscInt          i,n;
386c1a5e521SShri Abhyankar 
387c1a5e521SShri Abhyankar   PetscFunctionBegin;
388c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
389c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
390c1a5e521SShri Abhyankar   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
391c1a5e521SShri Abhyankar   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
392c1a5e521SShri Abhyankar 
393c1a5e521SShri Abhyankar   for(i = 0;i<n;i++) {
394c1a5e521SShri Abhyankar     if (x[i] < xl[i]) x[i] = xl[i];
395c1a5e521SShri Abhyankar     else if (x[i] > xu[i]) x[i] = xu[i];
396c1a5e521SShri Abhyankar   }
397c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
398c1a5e521SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
399c1a5e521SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
400c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
401c1a5e521SShri Abhyankar }
402c1a5e521SShri Abhyankar 
4032b4ed282SShri Abhyankar /*  --------------------------------------------------------------------
4042b4ed282SShri Abhyankar 
4052b4ed282SShri Abhyankar      This file implements a semismooth truncated Newton method with a line search,
4062b4ed282SShri Abhyankar      for solving a system of nonlinear equations in complementarity form, using the KSP, Vec,
4072b4ed282SShri Abhyankar      and Mat interfaces for linear solvers, vectors, and matrices,
4082b4ed282SShri Abhyankar      respectively.
4092b4ed282SShri Abhyankar 
4102b4ed282SShri Abhyankar      The following basic routines are required for each nonlinear solver:
4112b4ed282SShri Abhyankar           SNESCreate_XXX()          - Creates a nonlinear solver context
4122b4ed282SShri Abhyankar           SNESSetFromOptions_XXX()  - Sets runtime options
4132b4ed282SShri Abhyankar           SNESSolve_XXX()           - Solves the nonlinear system
4142b4ed282SShri Abhyankar           SNESDestroy_XXX()         - Destroys the nonlinear solver context
4152b4ed282SShri Abhyankar      The suffix "_XXX" denotes a particular implementation, in this case
4162b4ed282SShri Abhyankar      we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving
4172b4ed282SShri Abhyankar      systems of nonlinear equations with a line search (LS) method.
4182b4ed282SShri Abhyankar      These routines are actually called via the common user interface
4192b4ed282SShri Abhyankar      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
4202b4ed282SShri Abhyankar      SNESDestroy(), so the application code interface remains identical
4212b4ed282SShri Abhyankar      for all nonlinear solvers.
4222b4ed282SShri Abhyankar 
4232b4ed282SShri Abhyankar      Another key routine is:
4242b4ed282SShri Abhyankar           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
4252b4ed282SShri Abhyankar      by setting data structures and options.   The interface routine SNESSetUp()
4262b4ed282SShri Abhyankar      is not usually called directly by the user, but instead is called by
4272b4ed282SShri Abhyankar      SNESSolve() if necessary.
4282b4ed282SShri Abhyankar 
4292b4ed282SShri Abhyankar      Additional basic routines are:
4302b4ed282SShri Abhyankar           SNESView_XXX()            - Prints details of runtime options that
4312b4ed282SShri Abhyankar                                       have actually been used.
4322b4ed282SShri Abhyankar      These are called by application codes via the interface routines
4332b4ed282SShri Abhyankar      SNESView().
4342b4ed282SShri Abhyankar 
4352b4ed282SShri Abhyankar      The various types of solvers (preconditioners, Krylov subspace methods,
4362b4ed282SShri Abhyankar      nonlinear solvers, timesteppers) are all organized similarly, so the
4372b4ed282SShri Abhyankar      above description applies to these categories also.
4382b4ed282SShri Abhyankar 
4392b4ed282SShri Abhyankar     -------------------------------------------------------------------- */
4402b4ed282SShri Abhyankar /*
44169c03620SShri Abhyankar    SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton
4422b4ed282SShri Abhyankar    method using a line search.
4432b4ed282SShri Abhyankar 
4442b4ed282SShri Abhyankar    Input Parameters:
4452b4ed282SShri Abhyankar .  snes - the SNES context
4462b4ed282SShri Abhyankar 
4472b4ed282SShri Abhyankar    Output Parameter:
4482b4ed282SShri Abhyankar .  outits - number of iterations until termination
4492b4ed282SShri Abhyankar 
4502b4ed282SShri Abhyankar    Application Interface Routine: SNESSolve()
4512b4ed282SShri Abhyankar 
4522b4ed282SShri Abhyankar    Notes:
4532b4ed282SShri Abhyankar    This implements essentially a semismooth Newton method with a
45451defd01SShri Abhyankar    line search. The default line search does not do any line seach
45551defd01SShri Abhyankar    but rather takes a full newton step.
4562b4ed282SShri Abhyankar */
4572b4ed282SShri Abhyankar #undef __FUNCT__
45869c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS"
45969c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes)
4602b4ed282SShri Abhyankar {
4612b4ed282SShri Abhyankar   SNES_VI            *vi = (SNES_VI*)snes->data;
4622b4ed282SShri Abhyankar   PetscErrorCode     ierr;
4632b4ed282SShri Abhyankar   PetscInt           maxits,i,lits;
4643b336b49SShri Abhyankar   PetscBool          lssucceed;
4652b4ed282SShri Abhyankar   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
4662b4ed282SShri Abhyankar   PetscReal          gnorm,xnorm=0,ynorm;
4672b4ed282SShri Abhyankar   Vec                Y,X,F,G,W;
4682b4ed282SShri Abhyankar   KSPConvergedReason kspreason;
4692b4ed282SShri Abhyankar 
4702b4ed282SShri Abhyankar   PetscFunctionBegin;
4719ce7406fSBarry Smith   vi->computeuserfunction    = snes->ops->computefunction;
4729ce7406fSBarry Smith   snes->ops->computefunction = SNESVIComputeFunction;
4739ce7406fSBarry Smith 
4742b4ed282SShri Abhyankar   snes->numFailures            = 0;
4752b4ed282SShri Abhyankar   snes->numLinearSolveFailures = 0;
4762b4ed282SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
4772b4ed282SShri Abhyankar 
4782b4ed282SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
4792b4ed282SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
4802b4ed282SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
4812b4ed282SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
4822b4ed282SShri Abhyankar   G		= snes->work[1];
4832b4ed282SShri Abhyankar   W		= snes->work[2];
4842b4ed282SShri Abhyankar 
4852b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
4862b4ed282SShri Abhyankar   snes->iter = 0;
4872b4ed282SShri Abhyankar   snes->norm = 0.0;
4882b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
4892b4ed282SShri Abhyankar 
4907fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
4912b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
4922b4ed282SShri Abhyankar   if (snes->domainerror) {
4932b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
4949ce7406fSBarry Smith     snes->ops->computefunction = vi->computeuserfunction;
4952b4ed282SShri Abhyankar     PetscFunctionReturn(0);
4962b4ed282SShri Abhyankar   }
4972b4ed282SShri Abhyankar    /* Compute Merit function */
4982b4ed282SShri Abhyankar   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
4992b4ed282SShri Abhyankar 
5002b4ed282SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
5012b4ed282SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
5022b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5032b4ed282SShri Abhyankar 
5042b4ed282SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5052b4ed282SShri Abhyankar   snes->norm = vi->phinorm;
5062b4ed282SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5072b4ed282SShri Abhyankar   SNESLogConvHistory(snes,vi->phinorm,0);
5082b4ed282SShri Abhyankar   SNESMonitor(snes,0,vi->phinorm);
5092b4ed282SShri Abhyankar 
5102b4ed282SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
5112b4ed282SShri Abhyankar   snes->ttol = vi->phinorm*snes->rtol;
5122b4ed282SShri Abhyankar   /* test convergence */
5132b4ed282SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
5149ce7406fSBarry Smith   if (snes->reason) {
5159ce7406fSBarry Smith     snes->ops->computefunction = vi->computeuserfunction;
5169ce7406fSBarry Smith     PetscFunctionReturn(0);
5179ce7406fSBarry Smith   }
5182b4ed282SShri Abhyankar 
5192b4ed282SShri Abhyankar   for (i=0; i<maxits; i++) {
5202b4ed282SShri Abhyankar 
5212b4ed282SShri Abhyankar     /* Call general purpose update function */
5222b4ed282SShri Abhyankar     if (snes->ops->update) {
5232b4ed282SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
5242b4ed282SShri Abhyankar     }
5252b4ed282SShri Abhyankar 
5262b4ed282SShri Abhyankar     /* Solve J Y = Phi, where J is the semismooth jacobian */
527a79edbefSShri Abhyankar     /* Get the nonlinear function jacobian */
5282b4ed282SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
529a79edbefSShri Abhyankar     /* Get the diagonal shift and row scaling vectors */
530a79edbefSShri Abhyankar     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
531a79edbefSShri Abhyankar     /* Compute the semismooth jacobian */
532a79edbefSShri Abhyankar     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
533ee657d29SShri Abhyankar     /* Compute the merit function gradient */
534ee657d29SShri Abhyankar     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
5352b4ed282SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
5362b4ed282SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
5372b4ed282SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
5383b336b49SShri Abhyankar 
5393b336b49SShri Abhyankar     if (kspreason < 0) {
5402b4ed282SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
5412b4ed282SShri 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);
5422b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
5432b4ed282SShri Abhyankar         break;
5442b4ed282SShri Abhyankar       }
5452b4ed282SShri Abhyankar     }
5462b4ed282SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
5472b4ed282SShri Abhyankar     snes->linear_its += lits;
5482b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
5492b4ed282SShri Abhyankar     /*
5502b4ed282SShri Abhyankar     if (vi->precheckstep) {
551ace3abfcSBarry Smith       PetscBool changed_y = PETSC_FALSE;
5522b4ed282SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
5532b4ed282SShri Abhyankar     }
5542b4ed282SShri Abhyankar 
5552b4ed282SShri Abhyankar     if (PetscLogPrintInfo){
5562b4ed282SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
5572b4ed282SShri Abhyankar     }
5582b4ed282SShri Abhyankar     */
5592b4ed282SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
5602b4ed282SShri Abhyankar          Y <- X - lambda*Y
5612b4ed282SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
5622b4ed282SShri Abhyankar     */
5632b4ed282SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
5642b4ed282SShri Abhyankar     ynorm = 1; gnorm = vi->phinorm;
5652b4ed282SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
5662b4ed282SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
5672b4ed282SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
5682b4ed282SShri Abhyankar     if (snes->domainerror) {
5692b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
5709ce7406fSBarry Smith       snes->ops->computefunction = vi->computeuserfunction;
5712b4ed282SShri Abhyankar       PetscFunctionReturn(0);
5722b4ed282SShri Abhyankar     }
5732b4ed282SShri Abhyankar     if (!lssucceed) {
5742b4ed282SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
575ace3abfcSBarry Smith 	PetscBool ismin;
5762b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
5772b4ed282SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
5782b4ed282SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
5792b4ed282SShri Abhyankar         break;
5802b4ed282SShri Abhyankar       }
5812b4ed282SShri Abhyankar     }
5822b4ed282SShri Abhyankar     /* Update function and solution vectors */
5832b4ed282SShri Abhyankar     vi->phinorm = gnorm;
5842b4ed282SShri Abhyankar     vi->merit = 0.5*vi->phinorm*vi->phinorm;
5852b4ed282SShri Abhyankar     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
5862b4ed282SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
5872b4ed282SShri Abhyankar     /* Monitor convergence */
5882b4ed282SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5892b4ed282SShri Abhyankar     snes->iter = i+1;
5902b4ed282SShri Abhyankar     snes->norm = vi->phinorm;
5912b4ed282SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
5922b4ed282SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
5932b4ed282SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
5942b4ed282SShri Abhyankar     /* Test for convergence, xnorm = || X || */
5952b4ed282SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
5962b4ed282SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
5972b4ed282SShri Abhyankar     if (snes->reason) break;
5982b4ed282SShri Abhyankar   }
5992b4ed282SShri Abhyankar   if (i == maxits) {
6002b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
6012b4ed282SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
6022b4ed282SShri Abhyankar   }
6039ce7406fSBarry Smith   snes->ops->computefunction = vi->computeuserfunction;
6042b4ed282SShri Abhyankar   PetscFunctionReturn(0);
6052b4ed282SShri Abhyankar }
60669c03620SShri Abhyankar 
60769c03620SShri Abhyankar #undef __FUNCT__
608693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
609693ddf92SShri Abhyankar /*
610693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
611693ddf92SShri Abhyankar 
612693ddf92SShri Abhyankar    Input parameter
613693ddf92SShri Abhyankar .  snes - the SNES context
614693ddf92SShri Abhyankar .  X    - the snes solution vector
615693ddf92SShri Abhyankar .  F    - the nonlinear function vector
616693ddf92SShri Abhyankar 
617693ddf92SShri Abhyankar    Output parameter
618693ddf92SShri Abhyankar .  ISact - active set index set
619693ddf92SShri Abhyankar  */
620693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
621d950fb63SShri Abhyankar {
622d950fb63SShri Abhyankar   PetscErrorCode   ierr;
623693ddf92SShri Abhyankar   SNES_VI          *vi = (SNES_VI*)snes->data;
624693ddf92SShri Abhyankar   Vec               Xl=vi->xl,Xu=vi->xu;
625693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
626693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
627d950fb63SShri Abhyankar 
628d950fb63SShri Abhyankar   PetscFunctionBegin;
629d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
630d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
631fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
632fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
633fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
634fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
635693ddf92SShri Abhyankar   /* Compute active set size */
636d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
637693ddf92SShri 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_isact++;
638d950fb63SShri Abhyankar   }
639693ddf92SShri Abhyankar 
640d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
641d950fb63SShri Abhyankar 
642693ddf92SShri Abhyankar   /* Set active set indices */
643d950fb63SShri Abhyankar   for(i=0; i < nlocal; i++) {
644693ddf92SShri 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_act[i1++] = ilow+i;
645d950fb63SShri Abhyankar   }
646d950fb63SShri Abhyankar 
647693ddf92SShri Abhyankar    /* Create active set IS */
64862298a1eSBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
649d950fb63SShri Abhyankar 
650fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
651fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
652fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
653fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
654d950fb63SShri Abhyankar   PetscFunctionReturn(0);
655d950fb63SShri Abhyankar }
656d950fb63SShri Abhyankar 
657693ddf92SShri Abhyankar #undef __FUNCT__
658693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
659693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
660693ddf92SShri Abhyankar {
661693ddf92SShri Abhyankar   PetscErrorCode     ierr;
662693ddf92SShri Abhyankar 
663693ddf92SShri Abhyankar   PetscFunctionBegin;
664693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
665693ddf92SShri Abhyankar   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
666693ddf92SShri Abhyankar   PetscFunctionReturn(0);
667693ddf92SShri Abhyankar }
668693ddf92SShri Abhyankar 
669dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
670dbd914b8SShri Abhyankar #undef __FUNCT__
671fe835674SShri Abhyankar #define __FUNCT__ "SNESVICreateSubVectors"
672fe835674SShri Abhyankar PetscErrorCode SNESVICreateSubVectors(SNES snes,PetscInt n,Vec* newv)
673dbd914b8SShri Abhyankar {
674dbd914b8SShri Abhyankar   PetscErrorCode ierr;
675dbd914b8SShri Abhyankar   Vec            v;
676dbd914b8SShri Abhyankar 
677dbd914b8SShri Abhyankar   PetscFunctionBegin;
678dbd914b8SShri Abhyankar   ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr);
679dbd914b8SShri Abhyankar   ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr);
680dbd914b8SShri Abhyankar   ierr = VecSetFromOptions(v);CHKERRQ(ierr);
681dbd914b8SShri Abhyankar   *newv = v;
682dbd914b8SShri Abhyankar 
683dbd914b8SShri Abhyankar   PetscFunctionReturn(0);
684dbd914b8SShri Abhyankar }
685dbd914b8SShri Abhyankar 
686732bb160SShri Abhyankar /* Resets the snes PC and KSP when the active set sizes change */
687732bb160SShri Abhyankar #undef __FUNCT__
688732bb160SShri Abhyankar #define __FUNCT__ "SNESVIResetPCandKSP"
689732bb160SShri Abhyankar PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
690732bb160SShri Abhyankar {
691732bb160SShri Abhyankar   PetscErrorCode         ierr;
692732bb160SShri Abhyankar   KSP                    kspnew,snesksp;
693732bb160SShri Abhyankar   PC                     pcnew;
694732bb160SShri Abhyankar   const MatSolverPackage stype;
695dbd914b8SShri Abhyankar 
696732bb160SShri Abhyankar   PetscFunctionBegin;
697732bb160SShri Abhyankar   /* The active and inactive set sizes have changed so need to create a new snes->ksp object */
698732bb160SShri Abhyankar   ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
699732bb160SShri Abhyankar   ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr);
700732bb160SShri Abhyankar   /* Copy over snes->ksp info */
701732bb160SShri Abhyankar   kspnew->pc_side = snesksp->pc_side;
702732bb160SShri Abhyankar   kspnew->rtol    = snesksp->rtol;
703732bb160SShri Abhyankar   kspnew->abstol    = snesksp->abstol;
704732bb160SShri Abhyankar   kspnew->max_it  = snesksp->max_it;
705732bb160SShri Abhyankar   ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr);
706732bb160SShri Abhyankar   ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr);
707732bb160SShri Abhyankar   ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr);
708732bb160SShri Abhyankar   ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
709732bb160SShri Abhyankar   ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr);
710732bb160SShri Abhyankar   ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr);
711732bb160SShri Abhyankar   ierr = KSPDestroy(snesksp);CHKERRQ(ierr);
712732bb160SShri Abhyankar   snes->ksp = kspnew;
713732bb160SShri Abhyankar   ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr);
714732bb160SShri Abhyankar   ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);
715732bb160SShri Abhyankar   PetscFunctionReturn(0);
716732bb160SShri Abhyankar }
71769c03620SShri Abhyankar 
71869c03620SShri Abhyankar 
719fe835674SShri Abhyankar #undef __FUNCT__
720fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
721fe835674SShri Abhyankar PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscScalar *fnorm)
722fe835674SShri Abhyankar {
723fe835674SShri Abhyankar   PetscErrorCode    ierr;
724fe835674SShri Abhyankar   SNES_VI           *vi = (SNES_VI*)snes->data;
725fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
726fe835674SShri Abhyankar   PetscInt          i,n;
727fe835674SShri Abhyankar   PetscReal         rnorm;
728fe835674SShri Abhyankar 
729fe835674SShri Abhyankar   PetscFunctionBegin;
730fe835674SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
731fe835674SShri Abhyankar   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
732fe835674SShri Abhyankar   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
733fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
734fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
735fe835674SShri Abhyankar   rnorm = 0.0;
736fe835674SShri Abhyankar   for (i=0; i<n; i++) {
737fe835674SShri 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];
738fe835674SShri Abhyankar   }
739fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
740fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
741fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
742fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
743fe835674SShri Abhyankar   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
744fe835674SShri Abhyankar   *fnorm = sqrt(*fnorm);
745fe835674SShri Abhyankar   PetscFunctionReturn(0);
746fe835674SShri Abhyankar }
747fe835674SShri Abhyankar 
748fe835674SShri Abhyankar /* Variational Inequality solver using reduce space method. No semismooth algorithm is
749fe835674SShri Abhyankar    implemented in this algorithm. It basically identifies the active variables and does
750fe835674SShri Abhyankar    a linear solve on the inactive variables. */
751d950fb63SShri Abhyankar #undef __FUNCT__
752d950fb63SShri Abhyankar #define __FUNCT__ "SNESSolveVI_RS"
753d950fb63SShri Abhyankar PetscErrorCode SNESSolveVI_RS(SNES snes)
754d950fb63SShri Abhyankar {
755d950fb63SShri Abhyankar   SNES_VI          *vi = (SNES_VI*)snes->data;
756d950fb63SShri Abhyankar   PetscErrorCode    ierr;
757abcba341SShri Abhyankar   PetscInt          maxits,i,lits;
758d950fb63SShri Abhyankar   PetscBool         lssucceed;
759d950fb63SShri Abhyankar   MatStructure      flg = DIFFERENT_NONZERO_PATTERN;
760fe835674SShri Abhyankar   PetscReal         fnorm,gnorm,xnorm=0,ynorm;
761d950fb63SShri Abhyankar   Vec                Y,X,F,G,W;
762d950fb63SShri Abhyankar   KSPConvergedReason kspreason;
763d950fb63SShri Abhyankar 
764d950fb63SShri Abhyankar   PetscFunctionBegin;
765d950fb63SShri Abhyankar   snes->numFailures            = 0;
766d950fb63SShri Abhyankar   snes->numLinearSolveFailures = 0;
767d950fb63SShri Abhyankar   snes->reason                 = SNES_CONVERGED_ITERATING;
768d950fb63SShri Abhyankar 
769d950fb63SShri Abhyankar   maxits	= snes->max_its;	/* maximum number of iterations */
770d950fb63SShri Abhyankar   X		= snes->vec_sol;	/* solution vector */
771d950fb63SShri Abhyankar   F		= snes->vec_func;	/* residual vector */
772d950fb63SShri Abhyankar   Y		= snes->work[0];	/* work vectors */
773d950fb63SShri Abhyankar   G		= snes->work[1];
774d950fb63SShri Abhyankar   W		= snes->work[2];
775d950fb63SShri Abhyankar 
776d950fb63SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
777d950fb63SShri Abhyankar   snes->iter = 0;
778d950fb63SShri Abhyankar   snes->norm = 0.0;
779d950fb63SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
780d950fb63SShri Abhyankar 
7817fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
782fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
783d950fb63SShri Abhyankar   if (snes->domainerror) {
784d950fb63SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
785d950fb63SShri Abhyankar     PetscFunctionReturn(0);
786d950fb63SShri Abhyankar   }
787fe835674SShri Abhyankar   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
788d950fb63SShri Abhyankar   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
789d950fb63SShri Abhyankar   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
790fe835674SShri Abhyankar   if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
791d950fb63SShri Abhyankar 
792d950fb63SShri Abhyankar   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
793fe835674SShri Abhyankar   snes->norm = fnorm;
794d950fb63SShri Abhyankar   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
795fe835674SShri Abhyankar   SNESLogConvHistory(snes,fnorm,0);
796fe835674SShri Abhyankar   SNESMonitor(snes,0,fnorm);
797d950fb63SShri Abhyankar 
798d950fb63SShri Abhyankar   /* set parameter for default relative tolerance convergence test */
799fe835674SShri Abhyankar   snes->ttol = fnorm*snes->rtol;
800d950fb63SShri Abhyankar   /* test convergence */
801fe835674SShri Abhyankar   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
802d950fb63SShri Abhyankar   if (snes->reason) PetscFunctionReturn(0);
803d950fb63SShri Abhyankar 
804d950fb63SShri Abhyankar   for (i=0; i<maxits; i++) {
805d950fb63SShri Abhyankar 
806d950fb63SShri Abhyankar     IS         IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
807d950fb63SShri Abhyankar     VecScatter scat_act,scat_inact;
808abcba341SShri Abhyankar     PetscInt   nis_act,nis_inact;
809fe835674SShri Abhyankar     Vec        Y_act,Y_inact,F_inact;
810d950fb63SShri Abhyankar     Mat        jac_inact_inact,prejac_inact_inact;
81162298a1eSBarry Smith     IS         keptrows;
812abcba341SShri Abhyankar     PetscBool  isequal;
813d950fb63SShri Abhyankar 
814d950fb63SShri Abhyankar     /* Call general purpose update function */
815d950fb63SShri Abhyankar     if (snes->ops->update) {
816d950fb63SShri Abhyankar       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
817d950fb63SShri Abhyankar     }
818d950fb63SShri Abhyankar     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
81962298a1eSBarry Smith 
820d950fb63SShri Abhyankar     /* Create active and inactive index sets */
821693ddf92SShri Abhyankar     ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr);
822d950fb63SShri Abhyankar 
823abcba341SShri Abhyankar     /* Create inactive set submatrix */
82462298a1eSBarry Smith     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
825b3a44c85SBarry Smith     ierr = MatFindNonzeroRows(jac_inact_inact,&keptrows);CHKERRQ(ierr);
82662298a1eSBarry Smith     if (keptrows) {
82762298a1eSBarry Smith       PetscInt       cnt,*nrows,k;
82862298a1eSBarry Smith       const PetscInt *krows,*inact;
82927d4218bSShri Abhyankar       PetscInt       rstart=jac_inact_inact->rmap->rstart;
83062298a1eSBarry Smith 
83162298a1eSBarry Smith       ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
83262298a1eSBarry Smith       ierr = ISDestroy(IS_act);CHKERRQ(ierr);
83362298a1eSBarry Smith 
83462298a1eSBarry Smith       ierr = ISGetLocalSize(keptrows,&cnt);CHKERRQ(ierr);
83562298a1eSBarry Smith       ierr = ISGetIndices(keptrows,&krows);CHKERRQ(ierr);
83662298a1eSBarry Smith       ierr = ISGetIndices(IS_inact,&inact);CHKERRQ(ierr);
83762298a1eSBarry Smith       ierr = PetscMalloc(cnt*sizeof(PetscInt),&nrows);CHKERRQ(ierr);
83862298a1eSBarry Smith       for (k=0; k<cnt; k++) {
83927d4218bSShri Abhyankar         nrows[k] = inact[krows[k]-rstart];
84062298a1eSBarry Smith       }
84162298a1eSBarry Smith       ierr = ISRestoreIndices(keptrows,&krows);CHKERRQ(ierr);
84262298a1eSBarry Smith       ierr = ISRestoreIndices(IS_inact,&inact);CHKERRQ(ierr);
84362298a1eSBarry Smith       ierr = ISDestroy(keptrows);CHKERRQ(ierr);
84462298a1eSBarry Smith       ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
84562298a1eSBarry Smith 
84627d4218bSShri Abhyankar       ierr = ISCreateGeneral(PETSC_COMM_WORLD,cnt,nrows,PETSC_OWN_POINTER,&IS_inact);CHKERRQ(ierr);
84727d4218bSShri Abhyankar       ierr = ISComplement(IS_inact,F->map->rstart,F->map->rend,&IS_act);CHKERRQ(ierr);
84862298a1eSBarry Smith       ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
84962298a1eSBarry Smith     }
85062298a1eSBarry Smith 
851d950fb63SShri Abhyankar     /* Get sizes of active and inactive sets */
852d950fb63SShri Abhyankar     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
853d950fb63SShri Abhyankar     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
854d950fb63SShri Abhyankar 
855d950fb63SShri Abhyankar     /* Create active and inactive set vectors */
856fe835674SShri Abhyankar     ierr = SNESVICreateSubVectors(snes,nis_inact,&F_inact);CHKERRQ(ierr);
857fe835674SShri Abhyankar     ierr = SNESVICreateSubVectors(snes,nis_act,&Y_act);CHKERRQ(ierr);
858fe835674SShri Abhyankar     ierr = SNESVICreateSubVectors(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
859d950fb63SShri Abhyankar 
860d950fb63SShri Abhyankar     /* Create scatter contexts */
861d950fb63SShri Abhyankar     ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
862d950fb63SShri Abhyankar     ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
863d950fb63SShri Abhyankar 
864d950fb63SShri Abhyankar     /* Do a vec scatter to active and inactive set vectors */
865fe835674SShri Abhyankar     ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
866fe835674SShri Abhyankar     ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
867d950fb63SShri Abhyankar 
868d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
869d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
870d950fb63SShri Abhyankar 
871d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
872d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873d950fb63SShri Abhyankar 
874d950fb63SShri Abhyankar     /* Active set direction = 0 */
875d950fb63SShri Abhyankar     ierr = VecSet(Y_act,0);CHKERRQ(ierr);
876d950fb63SShri Abhyankar     if (snes->jacobian != snes->jacobian_pre) {
877d950fb63SShri Abhyankar       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
878d950fb63SShri Abhyankar     } else prejac_inact_inact = jac_inact_inact;
879d950fb63SShri Abhyankar 
880abcba341SShri Abhyankar     ierr = ISEqual(vi->IS_inact_prev,IS_inact,&isequal);CHKERRQ(ierr);
881abcba341SShri Abhyankar     if (!isequal) {
882732bb160SShri Abhyankar       ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr);
883d950fb63SShri Abhyankar     }
884d950fb63SShri Abhyankar 
885d950fb63SShri Abhyankar     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
886fe835674SShri Abhyankar     ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr);
887d950fb63SShri Abhyankar     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
888d950fb63SShri Abhyankar     if (kspreason < 0) {
889d950fb63SShri Abhyankar       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
890d950fb63SShri 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);
891d950fb63SShri Abhyankar         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
892d950fb63SShri Abhyankar         break;
893d950fb63SShri Abhyankar       }
894d950fb63SShri Abhyankar      }
895d950fb63SShri Abhyankar 
896d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
897d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
898d950fb63SShri Abhyankar     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
899d950fb63SShri Abhyankar     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
900d950fb63SShri Abhyankar 
901fe835674SShri Abhyankar     ierr = VecDestroy(F_inact);CHKERRQ(ierr);
902d950fb63SShri Abhyankar     ierr = VecDestroy(Y_act);CHKERRQ(ierr);
903d950fb63SShri Abhyankar     ierr = VecDestroy(Y_inact);CHKERRQ(ierr);
904d950fb63SShri Abhyankar     ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr);
905d950fb63SShri Abhyankar     ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr);
906d950fb63SShri Abhyankar     ierr = ISDestroy(IS_act);CHKERRQ(ierr);
907abcba341SShri Abhyankar     if (!isequal) {
908abcba341SShri Abhyankar       ierr = ISDestroy(vi->IS_inact_prev);CHKERRQ(ierr);
909abcba341SShri Abhyankar       ierr = ISDuplicate(IS_inact,&vi->IS_inact_prev);CHKERRQ(ierr);
910abcba341SShri Abhyankar     }
911d950fb63SShri Abhyankar     ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
912d950fb63SShri Abhyankar     ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
913d950fb63SShri Abhyankar     if (snes->jacobian != snes->jacobian_pre) {
914d950fb63SShri Abhyankar       ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr);
915d950fb63SShri Abhyankar     }
916d950fb63SShri Abhyankar     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
917d950fb63SShri Abhyankar     snes->linear_its += lits;
918d950fb63SShri Abhyankar     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
919d950fb63SShri Abhyankar     /*
920d950fb63SShri Abhyankar     if (vi->precheckstep) {
921d950fb63SShri Abhyankar       PetscBool changed_y = PETSC_FALSE;
922d950fb63SShri Abhyankar       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
923d950fb63SShri Abhyankar     }
924d950fb63SShri Abhyankar 
925d950fb63SShri Abhyankar     if (PetscLogPrintInfo){
926d950fb63SShri Abhyankar       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
927d950fb63SShri Abhyankar     }
928d950fb63SShri Abhyankar     */
929d950fb63SShri Abhyankar     /* Compute a (scaled) negative update in the line search routine:
930d950fb63SShri Abhyankar          Y <- X - lambda*Y
931d950fb63SShri Abhyankar        and evaluate G = function(Y) (depends on the line search).
932d950fb63SShri Abhyankar     */
933d950fb63SShri Abhyankar     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
934fe835674SShri Abhyankar     ynorm = 1; gnorm = fnorm;
935fe835674SShri Abhyankar     ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
936fe835674SShri Abhyankar     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
9372b4ed282SShri Abhyankar     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
9382b4ed282SShri Abhyankar     if (snes->domainerror) {
9392b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
9402b4ed282SShri Abhyankar       PetscFunctionReturn(0);
9412b4ed282SShri Abhyankar     }
9422b4ed282SShri Abhyankar     if (!lssucceed) {
9432b4ed282SShri Abhyankar       if (++snes->numFailures >= snes->maxFailures) {
9442b4ed282SShri Abhyankar 	PetscBool ismin;
9452b4ed282SShri Abhyankar         snes->reason = SNES_DIVERGED_LINE_SEARCH;
9462b4ed282SShri Abhyankar         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
9472b4ed282SShri Abhyankar         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
9482b4ed282SShri Abhyankar         break;
9492b4ed282SShri Abhyankar       }
9502b4ed282SShri Abhyankar     }
9512b4ed282SShri Abhyankar     /* Update function and solution vectors */
952fe835674SShri Abhyankar     fnorm = gnorm;
953fe835674SShri Abhyankar     ierr = VecCopy(G,F);CHKERRQ(ierr);
9542b4ed282SShri Abhyankar     ierr = VecCopy(W,X);CHKERRQ(ierr);
9552b4ed282SShri Abhyankar     /* Monitor convergence */
9562b4ed282SShri Abhyankar     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
9572b4ed282SShri Abhyankar     snes->iter = i+1;
958fe835674SShri Abhyankar     snes->norm = fnorm;
9592b4ed282SShri Abhyankar     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
9602b4ed282SShri Abhyankar     SNESLogConvHistory(snes,snes->norm,lits);
9612b4ed282SShri Abhyankar     SNESMonitor(snes,snes->iter,snes->norm);
9622b4ed282SShri Abhyankar     /* Test for convergence, xnorm = || X || */
9632b4ed282SShri Abhyankar     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
964fe835674SShri Abhyankar     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
9652b4ed282SShri Abhyankar     if (snes->reason) break;
9662b4ed282SShri Abhyankar   }
9672b4ed282SShri Abhyankar   if (i == maxits) {
9682b4ed282SShri Abhyankar     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
9692b4ed282SShri Abhyankar     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
9702b4ed282SShri Abhyankar   }
9712b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9722b4ed282SShri Abhyankar }
9732b4ed282SShri Abhyankar 
9742b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
9752b4ed282SShri Abhyankar /*
9762b4ed282SShri Abhyankar    SNESSetUp_VI - Sets up the internal data structures for the later use
9772b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
9782b4ed282SShri Abhyankar 
9792b4ed282SShri Abhyankar    Input Parameter:
9802b4ed282SShri Abhyankar .  snes - the SNES context
9812b4ed282SShri Abhyankar .  x - the solution vector
9822b4ed282SShri Abhyankar 
9832b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
9842b4ed282SShri Abhyankar 
9852b4ed282SShri Abhyankar    Notes:
9862b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
9872b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
9882b4ed282SShri Abhyankar    the call to SNESSolve().
9892b4ed282SShri Abhyankar  */
9902b4ed282SShri Abhyankar #undef __FUNCT__
9912b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
9922b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
9932b4ed282SShri Abhyankar {
9942b4ed282SShri Abhyankar   PetscErrorCode ierr;
9952b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*) snes->data;
9962b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
9972b4ed282SShri Abhyankar 
9982b4ed282SShri Abhyankar   PetscFunctionBegin;
9992b4ed282SShri Abhyankar   if (!snes->vec_sol_update) {
10002b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
10012b4ed282SShri Abhyankar     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
10022b4ed282SShri Abhyankar   }
10032b4ed282SShri Abhyankar   if (!snes->work) {
10042b4ed282SShri Abhyankar     snes->nwork = 3;
10052b4ed282SShri Abhyankar     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
10062b4ed282SShri Abhyankar     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
10072b4ed282SShri Abhyankar   }
10082b4ed282SShri Abhyankar 
10092b4ed282SShri Abhyankar   /* If the lower and upper bound on variables are not set, set it to
10102b4ed282SShri Abhyankar      -Inf and Inf */
10112b4ed282SShri Abhyankar   if (!vi->xl && !vi->xu) {
10122b4ed282SShri Abhyankar     vi->usersetxbounds = PETSC_FALSE;
10132b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr);
10142b4ed282SShri Abhyankar     ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr);
10152b4ed282SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr);
10162b4ed282SShri Abhyankar     ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr);
10172b4ed282SShri Abhyankar   } else {
10182b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
10192b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
10202b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
10212b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
10222b4ed282SShri 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]))
10232b4ed282SShri 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.");
10242b4ed282SShri Abhyankar   }
1025abcba341SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
1026abcba341SShri Abhyankar     /* Set up previous active index set for the first snes solve
1027abcba341SShri Abhyankar        vi->IS_inact_prev = 0,1,2,....N */
1028abcba341SShri Abhyankar     PetscInt *indices;
1029abcba341SShri Abhyankar     PetscInt  i,n,rstart,rend;
1030abcba341SShri Abhyankar 
1031abcba341SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);CHKERRQ(ierr);
1032abcba341SShri Abhyankar     ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1033abcba341SShri Abhyankar     ierr = PetscMalloc(n*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1034abcba341SShri Abhyankar     for(i=0;i < n; i++) indices[i] = rstart + i;
1035abcba341SShri Abhyankar     ierr = ISCreateGeneral(((PetscObject)snes)->comm,n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);
1036abcba341SShri Abhyankar   }
10372b4ed282SShri Abhyankar 
1038fe835674SShri Abhyankar   if (snes->ops->solve != SNESSolveVI_RS) {
1039fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
1040fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr);
1041fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr);
1042fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr);
1043fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
1044fe835674SShri Abhyankar     ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr);
1045fe835674SShri Abhyankar   }
10462b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10472b4ed282SShri Abhyankar }
10482b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
10492b4ed282SShri Abhyankar /*
10502b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
10512b4ed282SShri Abhyankar    with SNESCreate_VI().
10522b4ed282SShri Abhyankar 
10532b4ed282SShri Abhyankar    Input Parameter:
10542b4ed282SShri Abhyankar .  snes - the SNES context
10552b4ed282SShri Abhyankar 
10562b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
10572b4ed282SShri Abhyankar  */
10582b4ed282SShri Abhyankar #undef __FUNCT__
10592b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
10602b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
10612b4ed282SShri Abhyankar {
10622b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*) snes->data;
10632b4ed282SShri Abhyankar   PetscErrorCode ierr;
10642b4ed282SShri Abhyankar 
10652b4ed282SShri Abhyankar   PetscFunctionBegin;
10662b4ed282SShri Abhyankar   if (snes->vec_sol_update) {
10672b4ed282SShri Abhyankar     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
10682b4ed282SShri Abhyankar     snes->vec_sol_update = PETSC_NULL;
10692b4ed282SShri Abhyankar   }
10702b4ed282SShri Abhyankar   if (snes->nwork) {
1071*1130f37dSBarry Smith     ierr = VecDestroyVecs(&snes->work,snes->nwork);CHKERRQ(ierr);
10722b4ed282SShri Abhyankar     snes->nwork = 0;
10732b4ed282SShri Abhyankar   }
1074abcba341SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
1075abcba341SShri Abhyankar     ierr = ISDestroy(vi->IS_inact_prev);
1076abcba341SShri Abhyankar   }
10772b4ed282SShri Abhyankar 
1078fe835674SShri Abhyankar   if (snes->ops->solve != SNESSolveVI_RS) {
10792b4ed282SShri Abhyankar     /* clear vectors */
1080ee657d29SShri Abhyankar     ierr = VecDestroy(vi->dpsi);CHKERRQ(ierr);
10812b4ed282SShri Abhyankar     ierr = VecDestroy(vi->phi); CHKERRQ(ierr);
10822b4ed282SShri Abhyankar     ierr = VecDestroy(vi->Da); CHKERRQ(ierr);
10832b4ed282SShri Abhyankar     ierr = VecDestroy(vi->Db); CHKERRQ(ierr);
10842b4ed282SShri Abhyankar     ierr = VecDestroy(vi->z); CHKERRQ(ierr);
10852b4ed282SShri Abhyankar     ierr = VecDestroy(vi->t); CHKERRQ(ierr);
10867fe79bc4SShri Abhyankar   }
10877fe79bc4SShri Abhyankar 
10882b4ed282SShri Abhyankar   if (!vi->usersetxbounds) {
10892b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xl); CHKERRQ(ierr);
10902b4ed282SShri Abhyankar     ierr = VecDestroy(vi->xu); CHKERRQ(ierr);
10912b4ed282SShri Abhyankar   }
10927fe79bc4SShri Abhyankar 
1093c92abb8aSShri Abhyankar   if (vi->lsmonitor) {
1094c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
10952b4ed282SShri Abhyankar   }
10962b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
10972b4ed282SShri Abhyankar 
10982b4ed282SShri Abhyankar   /* clear composed functions */
10992b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
1100c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
11012b4ed282SShri Abhyankar   PetscFunctionReturn(0);
11022b4ed282SShri Abhyankar }
11037fe79bc4SShri Abhyankar 
11042b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
11052b4ed282SShri Abhyankar #undef __FUNCT__
11062b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
11072b4ed282SShri Abhyankar 
11082b4ed282SShri Abhyankar /*
11097fe79bc4SShri Abhyankar   This routine does not actually do a line search but it takes a full newton
11107fe79bc4SShri Abhyankar   step while ensuring that the new iterates remain within the constraints.
11112b4ed282SShri Abhyankar 
11122b4ed282SShri Abhyankar */
1113ace3abfcSBarry 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)
11142b4ed282SShri Abhyankar {
11152b4ed282SShri Abhyankar   PetscErrorCode ierr;
11162b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1117ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
11182b4ed282SShri Abhyankar 
11192b4ed282SShri Abhyankar   PetscFunctionBegin;
11202b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
11212b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11222b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
11232b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1124c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1125c1a5e521SShri Abhyankar   if (vi->postcheckstep) {
1126c1a5e521SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1127c1a5e521SShri Abhyankar   }
1128c1a5e521SShri Abhyankar   if (changed_y) {
1129c1a5e521SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
11307fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1131c1a5e521SShri Abhyankar   }
1132c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1133c1a5e521SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1134c1a5e521SShri Abhyankar   if (!snes->domainerror) {
11357fe79bc4SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_RS) {
11367fe79bc4SShri Abhyankar        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
11377fe79bc4SShri Abhyankar     } else {
1138c1a5e521SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
11397fe79bc4SShri Abhyankar     }
1140c1a5e521SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1141c1a5e521SShri Abhyankar   }
1142c1a5e521SShri Abhyankar   if (vi->lsmonitor) {
1143c1a5e521SShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1144c1a5e521SShri Abhyankar   }
1145c1a5e521SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1146c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
1147c1a5e521SShri Abhyankar }
1148c1a5e521SShri Abhyankar 
1149c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
1150c1a5e521SShri Abhyankar #undef __FUNCT__
11512b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
11522b4ed282SShri Abhyankar 
11532b4ed282SShri Abhyankar /*
11542b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
11552b4ed282SShri Abhyankar */
1156ace3abfcSBarry 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)
11572b4ed282SShri Abhyankar {
11582b4ed282SShri Abhyankar   PetscErrorCode ierr;
11592b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1160ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
11612b4ed282SShri Abhyankar 
11622b4ed282SShri Abhyankar   PetscFunctionBegin;
11632b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
11642b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11652b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
11667fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
11672b4ed282SShri Abhyankar   if (vi->postcheckstep) {
11682b4ed282SShri Abhyankar    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
11692b4ed282SShri Abhyankar   }
11702b4ed282SShri Abhyankar   if (changed_y) {
11712b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
11727fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
11732b4ed282SShri Abhyankar   }
11742b4ed282SShri Abhyankar 
11752b4ed282SShri Abhyankar   /* don't evaluate function the last time through */
11762b4ed282SShri Abhyankar   if (snes->iter < snes->max_its-1) {
11772b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
11782b4ed282SShri Abhyankar   }
11792b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
11802b4ed282SShri Abhyankar   PetscFunctionReturn(0);
11812b4ed282SShri Abhyankar }
11827fe79bc4SShri Abhyankar 
11832b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
11842b4ed282SShri Abhyankar #undef __FUNCT__
11852b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
11862b4ed282SShri Abhyankar /*
11877fe79bc4SShri Abhyankar   This routine implements a cubic line search while doing a projection on the variable bounds
11882b4ed282SShri Abhyankar */
1189ace3abfcSBarry 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)
11902b4ed282SShri Abhyankar {
1191fe835674SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
1192fe835674SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
1193fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
1194fe835674SShri Abhyankar   PetscScalar    cinitslope;
1195fe835674SShri Abhyankar #endif
1196fe835674SShri Abhyankar   PetscErrorCode ierr;
1197fe835674SShri Abhyankar   PetscInt       count;
1198fe835674SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1199fe835674SShri Abhyankar   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1200fe835674SShri Abhyankar   MPI_Comm       comm;
1201fe835674SShri Abhyankar 
1202fe835674SShri Abhyankar   PetscFunctionBegin;
1203fe835674SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1204fe835674SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1205fe835674SShri Abhyankar   *flag   = PETSC_TRUE;
1206fe835674SShri Abhyankar 
1207fe835674SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1208fe835674SShri Abhyankar   if (*ynorm == 0.0) {
1209fe835674SShri Abhyankar     if (vi->lsmonitor) {
1210fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
1211fe835674SShri Abhyankar     }
1212fe835674SShri Abhyankar     *gnorm = fnorm;
1213fe835674SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1214fe835674SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1215fe835674SShri Abhyankar     *flag  = PETSC_FALSE;
1216fe835674SShri Abhyankar     goto theend1;
1217fe835674SShri Abhyankar   }
1218fe835674SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1219fe835674SShri Abhyankar     if (vi->lsmonitor) {
1220fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
1221fe835674SShri Abhyankar     }
1222fe835674SShri Abhyankar     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1223fe835674SShri Abhyankar     *ynorm = vi->maxstep;
1224fe835674SShri Abhyankar   }
1225fe835674SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1226fe835674SShri Abhyankar   minlambda = vi->minlambda/rellength;
1227fe835674SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1228fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
1229fe835674SShri Abhyankar   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1230fe835674SShri Abhyankar   initslope = PetscRealPart(cinitslope);
1231fe835674SShri Abhyankar #else
1232fe835674SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1233fe835674SShri Abhyankar #endif
1234fe835674SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
1235fe835674SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
1236fe835674SShri Abhyankar 
1237fe835674SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1238fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1239fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
1240fe835674SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1241fe835674SShri Abhyankar     *flag = PETSC_FALSE;
1242fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1243fe835674SShri Abhyankar     goto theend1;
1244fe835674SShri Abhyankar   }
1245fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12467fe79bc4SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
1247fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
12487fe79bc4SShri Abhyankar   } else {
12497fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12507fe79bc4SShri Abhyankar   }
1251fe835674SShri Abhyankar   if (snes->domainerror) {
1252fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1253fe835674SShri Abhyankar     PetscFunctionReturn(0);
1254fe835674SShri Abhyankar   }
1255fe835674SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1256fe835674SShri Abhyankar   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1257fe835674SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1258fe835674SShri Abhyankar     if (vi->lsmonitor) {
1259fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1260fe835674SShri Abhyankar     }
1261fe835674SShri Abhyankar     goto theend1;
1262fe835674SShri Abhyankar   }
1263fe835674SShri Abhyankar 
1264fe835674SShri Abhyankar   /* Fit points with quadratic */
1265fe835674SShri Abhyankar   lambda     = 1.0;
1266fe835674SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
1267fe835674SShri Abhyankar   lambdaprev = lambda;
1268fe835674SShri Abhyankar   gnormprev  = *gnorm;
1269fe835674SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1270fe835674SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
1271fe835674SShri Abhyankar   else                         lambda = lambdatemp;
1272fe835674SShri Abhyankar 
1273fe835674SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1274fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1275fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
1276fe835674SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
1277fe835674SShri Abhyankar     *flag = PETSC_FALSE;
1278fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1279fe835674SShri Abhyankar     goto theend1;
1280fe835674SShri Abhyankar   }
1281fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
12827fe79bc4SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
1283fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
12847fe79bc4SShri Abhyankar   } else {
12857fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
12867fe79bc4SShri Abhyankar   }
1287fe835674SShri Abhyankar   if (snes->domainerror) {
1288fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1289fe835674SShri Abhyankar     PetscFunctionReturn(0);
1290fe835674SShri Abhyankar   }
1291fe835674SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1292fe835674SShri Abhyankar   if (vi->lsmonitor) {
1293fe835674SShri Abhyankar     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
1294fe835674SShri Abhyankar   }
1295fe835674SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1296fe835674SShri Abhyankar     if (vi->lsmonitor) {
1297fe835674SShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
1298fe835674SShri Abhyankar     }
1299fe835674SShri Abhyankar     goto theend1;
1300fe835674SShri Abhyankar   }
1301fe835674SShri Abhyankar 
1302fe835674SShri Abhyankar   /* Fit points with cubic */
1303fe835674SShri Abhyankar   count = 1;
1304fe835674SShri Abhyankar   while (PETSC_TRUE) {
1305fe835674SShri Abhyankar     if (lambda <= minlambda) {
1306fe835674SShri Abhyankar       if (vi->lsmonitor) {
1307fe835674SShri Abhyankar 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
1308fe835674SShri 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);
1309fe835674SShri Abhyankar       }
1310fe835674SShri Abhyankar       *flag = PETSC_FALSE;
1311fe835674SShri Abhyankar       break;
1312fe835674SShri Abhyankar     }
1313fe835674SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
1314fe835674SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
1315fe835674SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1316fe835674SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1317fe835674SShri Abhyankar     d  = b*b - 3*a*initslope;
1318fe835674SShri Abhyankar     if (d < 0.0) d = 0.0;
1319fe835674SShri Abhyankar     if (a == 0.0) {
1320fe835674SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
1321fe835674SShri Abhyankar     } else {
1322fe835674SShri Abhyankar       lambdatemp = (-b + sqrt(d))/(3.0*a);
1323fe835674SShri Abhyankar     }
1324fe835674SShri Abhyankar     lambdaprev = lambda;
1325fe835674SShri Abhyankar     gnormprev  = *gnorm;
1326fe835674SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1327fe835674SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
1328fe835674SShri Abhyankar     else                         lambda     = lambdatemp;
1329fe835674SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1330fe835674SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1331fe835674SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
1332fe835674SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
1333fe835674SShri 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);
1334fe835674SShri Abhyankar       *flag = PETSC_FALSE;
1335fe835674SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1336fe835674SShri Abhyankar       break;
1337fe835674SShri Abhyankar     }
1338fe835674SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13397fe79bc4SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_RS) {
1340fe835674SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
13417fe79bc4SShri Abhyankar     } else {
13427fe79bc4SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13437fe79bc4SShri Abhyankar     }
1344fe835674SShri Abhyankar     if (snes->domainerror) {
1345fe835674SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1346fe835674SShri Abhyankar       PetscFunctionReturn(0);
1347fe835674SShri Abhyankar     }
1348fe835674SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1349fe835674SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
1350fe835674SShri Abhyankar       if (vi->lsmonitor) {
1351fe835674SShri Abhyankar 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1352fe835674SShri Abhyankar       }
1353fe835674SShri Abhyankar       break;
1354fe835674SShri Abhyankar     } else {
1355fe835674SShri Abhyankar       if (vi->lsmonitor) {
1356fe835674SShri Abhyankar         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1357fe835674SShri Abhyankar       }
1358fe835674SShri Abhyankar     }
1359fe835674SShri Abhyankar     count++;
1360fe835674SShri Abhyankar   }
1361fe835674SShri Abhyankar   theend1:
1362fe835674SShri Abhyankar   /* Optional user-defined check for line search step validity */
1363fe835674SShri Abhyankar   if (vi->postcheckstep && *flag) {
1364fe835674SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1365fe835674SShri Abhyankar     if (changed_y) {
1366fe835674SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1367fe835674SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1368fe835674SShri Abhyankar     }
1369fe835674SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
1370fe835674SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
13717fe79bc4SShri Abhyankar       if (snes->ops->solve == SNESSolveVI_RS) {
1372fe835674SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
13737fe79bc4SShri Abhyankar       } else {
13747fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
13757fe79bc4SShri Abhyankar       }
1376fe835674SShri Abhyankar       if (snes->domainerror) {
1377fe835674SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1378fe835674SShri Abhyankar         PetscFunctionReturn(0);
1379fe835674SShri Abhyankar       }
1380fe835674SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1381fe835674SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
1382fe835674SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
1383fe835674SShri Abhyankar     }
1384fe835674SShri Abhyankar   }
1385fe835674SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1386fe835674SShri Abhyankar   PetscFunctionReturn(0);
1387fe835674SShri Abhyankar }
1388fe835674SShri Abhyankar 
13892b4ed282SShri Abhyankar #undef __FUNCT__
13902b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
13912b4ed282SShri Abhyankar /*
13927fe79bc4SShri Abhyankar   This routine does a quadratic line search while keeping the iterates within the variable bounds
13932b4ed282SShri Abhyankar */
1394ace3abfcSBarry 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)
13952b4ed282SShri Abhyankar {
13962b4ed282SShri Abhyankar   /*
13972b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
13982b4ed282SShri Abhyankar      minimization problem:
13992b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
14002b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
14012b4ed282SShri Abhyankar    */
14022b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
14032b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
14042b4ed282SShri Abhyankar   PetscScalar    cinitslope;
14052b4ed282SShri Abhyankar #endif
14062b4ed282SShri Abhyankar   PetscErrorCode ierr;
14072b4ed282SShri Abhyankar   PetscInt       count;
14082b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
1409ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
14102b4ed282SShri Abhyankar 
14112b4ed282SShri Abhyankar   PetscFunctionBegin;
14122b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14132b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
14142b4ed282SShri Abhyankar 
14152b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
14162b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
1417c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1418c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
14192b4ed282SShri Abhyankar     }
14202b4ed282SShri Abhyankar     *gnorm = fnorm;
14212b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
14222b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
14232b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
14242b4ed282SShri Abhyankar     goto theend2;
14252b4ed282SShri Abhyankar   }
14262b4ed282SShri Abhyankar   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
14272b4ed282SShri Abhyankar     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
14282b4ed282SShri Abhyankar     *ynorm = vi->maxstep;
14292b4ed282SShri Abhyankar   }
14302b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
14312b4ed282SShri Abhyankar   minlambda = vi->minlambda/rellength;
14327fe79bc4SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
14332b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
14347fe79bc4SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
14352b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
14362b4ed282SShri Abhyankar #else
14377fe79bc4SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
14382b4ed282SShri Abhyankar #endif
14392b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
14402b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
14412b4ed282SShri Abhyankar   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
14422b4ed282SShri Abhyankar 
14432b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
14447fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
14452b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
14462b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
14472b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
14482b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
14492b4ed282SShri Abhyankar     goto theend2;
14502b4ed282SShri Abhyankar   }
14512b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
14527fe79bc4SShri Abhyankar   if (snes->ops->solve == SNESSolveVI_RS) {
14537fe79bc4SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
14547fe79bc4SShri Abhyankar   } else {
14557fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
14567fe79bc4SShri Abhyankar   }
14572b4ed282SShri Abhyankar   if (snes->domainerror) {
14582b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14592b4ed282SShri Abhyankar     PetscFunctionReturn(0);
14602b4ed282SShri Abhyankar   }
14612b4ed282SShri Abhyankar   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
14622b4ed282SShri Abhyankar   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1463c92abb8aSShri Abhyankar     if (vi->lsmonitor) {
1464c92abb8aSShri Abhyankar       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
14652b4ed282SShri Abhyankar     }
14662b4ed282SShri Abhyankar     goto theend2;
14672b4ed282SShri Abhyankar   }
14682b4ed282SShri Abhyankar 
14692b4ed282SShri Abhyankar   /* Fit points with quadratic */
14702b4ed282SShri Abhyankar   lambda = 1.0;
14712b4ed282SShri Abhyankar   count = 1;
14722b4ed282SShri Abhyankar   while (PETSC_TRUE) {
14732b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
1474c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1475c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
1476c92abb8aSShri 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);
14772b4ed282SShri Abhyankar       }
14782b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
14792b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
14802b4ed282SShri Abhyankar       break;
14812b4ed282SShri Abhyankar     }
14822b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
14832b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
14842b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
14852b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
14862b4ed282SShri Abhyankar 
14872b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
14887fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
14892b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
14902b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
14912b4ed282SShri 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);
14922b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
14932b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
14942b4ed282SShri Abhyankar       break;
14952b4ed282SShri Abhyankar     }
14962b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
14972b4ed282SShri Abhyankar     if (snes->domainerror) {
14982b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
14992b4ed282SShri Abhyankar       PetscFunctionReturn(0);
15002b4ed282SShri Abhyankar     }
15017fe79bc4SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_RS) {
15027fe79bc4SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
15037fe79bc4SShri Abhyankar     } else {
15042b4ed282SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
15057fe79bc4SShri Abhyankar     }
15062b4ed282SShri Abhyankar     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
15072b4ed282SShri Abhyankar     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1508c92abb8aSShri Abhyankar       if (vi->lsmonitor) {
1509c92abb8aSShri Abhyankar         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
15102b4ed282SShri Abhyankar       }
15112b4ed282SShri Abhyankar       break;
15122b4ed282SShri Abhyankar     }
15132b4ed282SShri Abhyankar     count++;
15142b4ed282SShri Abhyankar   }
15152b4ed282SShri Abhyankar   theend2:
15162b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
15172b4ed282SShri Abhyankar   if (vi->postcheckstep) {
15182b4ed282SShri Abhyankar     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
15192b4ed282SShri Abhyankar     if (changed_y) {
15202b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
15217fe79bc4SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
15222b4ed282SShri Abhyankar     }
15232b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
15242b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
15252b4ed282SShri Abhyankar       if (snes->domainerror) {
15262b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
15272b4ed282SShri Abhyankar         PetscFunctionReturn(0);
15282b4ed282SShri Abhyankar       }
15297fe79bc4SShri Abhyankar       if (snes->ops->solve == SNESSolveVI_RS) {
15307fe79bc4SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
15317fe79bc4SShri Abhyankar       } else {
15327fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
15337fe79bc4SShri Abhyankar       }
15347fe79bc4SShri Abhyankar 
15352b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
15362b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
15372b4ed282SShri Abhyankar       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
15382b4ed282SShri Abhyankar     }
15392b4ed282SShri Abhyankar   }
15402b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
15412b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15422b4ed282SShri Abhyankar }
15432b4ed282SShri Abhyankar 
1544ace3abfcSBarry 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*/
15452b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15462b4ed282SShri Abhyankar EXTERN_C_BEGIN
15472b4ed282SShri Abhyankar #undef __FUNCT__
15482b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI"
15497087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
15502b4ed282SShri Abhyankar {
15512b4ed282SShri Abhyankar   PetscFunctionBegin;
15522b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->LineSearch = func;
15532b4ed282SShri Abhyankar   ((SNES_VI *)(snes->data))->lsP        = lsctx;
15542b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15552b4ed282SShri Abhyankar }
15562b4ed282SShri Abhyankar EXTERN_C_END
15572b4ed282SShri Abhyankar 
15582b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
15592b4ed282SShri Abhyankar EXTERN_C_BEGIN
15602b4ed282SShri Abhyankar #undef __FUNCT__
15612b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
15627087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
15632b4ed282SShri Abhyankar {
15642b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI*)snes->data;
15652b4ed282SShri Abhyankar   PetscErrorCode ierr;
15662b4ed282SShri Abhyankar 
15672b4ed282SShri Abhyankar   PetscFunctionBegin;
1568c92abb8aSShri Abhyankar   if (flg && !vi->lsmonitor) {
1569c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr);
1570c92abb8aSShri Abhyankar   } else if (!flg && vi->lsmonitor) {
1571c92abb8aSShri Abhyankar     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
15722b4ed282SShri Abhyankar   }
15732b4ed282SShri Abhyankar   PetscFunctionReturn(0);
15742b4ed282SShri Abhyankar }
15752b4ed282SShri Abhyankar EXTERN_C_END
15762b4ed282SShri Abhyankar 
15772b4ed282SShri Abhyankar /*
15782b4ed282SShri Abhyankar    SNESView_VI - Prints info from the SNESVI data structure.
15792b4ed282SShri Abhyankar 
15802b4ed282SShri Abhyankar    Input Parameters:
15812b4ed282SShri Abhyankar .  SNES - the SNES context
15822b4ed282SShri Abhyankar .  viewer - visualization context
15832b4ed282SShri Abhyankar 
15842b4ed282SShri Abhyankar    Application Interface Routine: SNESView()
15852b4ed282SShri Abhyankar */
15862b4ed282SShri Abhyankar #undef __FUNCT__
15872b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI"
15882b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
15892b4ed282SShri Abhyankar {
15902b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
159178c4b9d3SShri Abhyankar   const char     *cstr,*tstr;
15922b4ed282SShri Abhyankar   PetscErrorCode ierr;
1593ace3abfcSBarry Smith   PetscBool     iascii;
15942b4ed282SShri Abhyankar 
15952b4ed282SShri Abhyankar   PetscFunctionBegin;
15962b4ed282SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
15972b4ed282SShri Abhyankar   if (iascii) {
15982b4ed282SShri Abhyankar     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
15992b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
16002b4ed282SShri Abhyankar     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
16012b4ed282SShri Abhyankar     else                                                             cstr = "unknown";
160278c4b9d3SShri Abhyankar     if (snes->ops->solve == SNESSolveVI_SS)      tstr = "Semismooth";
16030a7c7af0SShri Abhyankar     else if (snes->ops->solve == SNESSolveVI_RS) tstr = "Reduced Space";
160478c4b9d3SShri Abhyankar     else                                         tstr = "unknown";
160578c4b9d3SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  VI algorithm: %s\n",tstr);CHKERRQ(ierr);
16062b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
16072b4ed282SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
16082b4ed282SShri Abhyankar   } else {
16092b4ed282SShri Abhyankar     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
16102b4ed282SShri Abhyankar   }
16112b4ed282SShri Abhyankar   PetscFunctionReturn(0);
16122b4ed282SShri Abhyankar }
16132b4ed282SShri Abhyankar 
16145559b345SBarry Smith #undef __FUNCT__
16155559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
16165559b345SBarry Smith /*@
16172b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
16182b4ed282SShri Abhyankar 
16192b4ed282SShri Abhyankar    Input Parameters:
16202b4ed282SShri Abhyankar .  snes - the SNES context.
16212b4ed282SShri Abhyankar .  xl   - lower bound.
16222b4ed282SShri Abhyankar .  xu   - upper bound.
16232b4ed282SShri Abhyankar 
16242b4ed282SShri Abhyankar    Notes:
16252b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
162684c105d7SBarry Smith    PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp().
162784c105d7SBarry Smith 
16285559b345SBarry Smith @*/
162969c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
16302b4ed282SShri Abhyankar {
16312b4ed282SShri Abhyankar   SNES_VI  *vi = (SNES_VI*)snes->data;
16322b4ed282SShri Abhyankar 
16332b4ed282SShri Abhyankar   PetscFunctionBegin;
16342b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
16352b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
16362b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
16372b4ed282SShri Abhyankar   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
16382b4ed282SShri 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);
16392b4ed282SShri 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);
16405559b345SBarry Smith 
16412b4ed282SShri Abhyankar   vi->usersetxbounds = PETSC_TRUE;
16422b4ed282SShri Abhyankar   vi->xl = xl;
16432b4ed282SShri Abhyankar   vi->xu = xu;
16442b4ed282SShri Abhyankar   PetscFunctionReturn(0);
16452b4ed282SShri Abhyankar }
1646693ddf92SShri Abhyankar 
16472b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
16482b4ed282SShri Abhyankar /*
16492b4ed282SShri Abhyankar    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
16502b4ed282SShri Abhyankar 
16512b4ed282SShri Abhyankar    Input Parameter:
16522b4ed282SShri Abhyankar .  snes - the SNES context
16532b4ed282SShri Abhyankar 
16542b4ed282SShri Abhyankar    Application Interface Routine: SNESSetFromOptions()
16552b4ed282SShri Abhyankar */
16562b4ed282SShri Abhyankar #undef __FUNCT__
16572b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI"
16582b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
16592b4ed282SShri Abhyankar {
16602b4ed282SShri Abhyankar   SNES_VI        *vi = (SNES_VI *)snes->data;
16617fe79bc4SShri Abhyankar   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
16625559b345SBarry Smith   const char     *vies[] = {"ss","rs"};
16632b4ed282SShri Abhyankar   PetscErrorCode ierr;
16642b4ed282SShri Abhyankar   PetscInt       indx;
166569c03620SShri Abhyankar   PetscBool      flg,set,flg2;
16662b4ed282SShri Abhyankar 
16672b4ed282SShri Abhyankar   PetscFunctionBegin;
16682b4ed282SShri Abhyankar   ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
16699308c137SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
16709308c137SBarry Smith   if (flg) {
16719308c137SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
16729308c137SBarry Smith   }
1673be6adb11SBarry Smith   ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
1674be6adb11SBarry Smith   ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
1675be6adb11SBarry Smith   ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
16762b4ed282SShri Abhyankar   ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
1677be6adb11SBarry Smith   ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
16782b4ed282SShri Abhyankar   if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
16795559b345SBarry Smith   ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr);
168069c03620SShri Abhyankar   if (flg2) {
168169c03620SShri Abhyankar     switch (indx) {
168269c03620SShri Abhyankar     case 0:
168369c03620SShri Abhyankar       snes->ops->solve = SNESSolveVI_SS;
168469c03620SShri Abhyankar       break;
168569c03620SShri Abhyankar     case 1:
1686d950fb63SShri Abhyankar       snes->ops->solve = SNESSolveVI_RS;
1687d950fb63SShri Abhyankar       break;
168869c03620SShri Abhyankar     }
168969c03620SShri Abhyankar   }
1690be6adb11SBarry Smith   ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"basic",&indx,&flg);CHKERRQ(ierr);
16912b4ed282SShri Abhyankar   if (flg) {
16922b4ed282SShri Abhyankar     switch (indx) {
16932b4ed282SShri Abhyankar     case 0:
16942b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
16952b4ed282SShri Abhyankar       break;
16962b4ed282SShri Abhyankar     case 1:
16972b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
16982b4ed282SShri Abhyankar       break;
16992b4ed282SShri Abhyankar     case 2:
17002b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
17012b4ed282SShri Abhyankar       break;
17022b4ed282SShri Abhyankar     case 3:
17032b4ed282SShri Abhyankar       ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
17042b4ed282SShri Abhyankar       break;
17052b4ed282SShri Abhyankar     }
1706fe835674SShri Abhyankar   }
17072b4ed282SShri Abhyankar   ierr = PetscOptionsTail();CHKERRQ(ierr);
17082b4ed282SShri Abhyankar   PetscFunctionReturn(0);
17092b4ed282SShri Abhyankar }
17102b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
17112b4ed282SShri Abhyankar /*MC
17122b4ed282SShri Abhyankar       SNESVI - Semismooth newton method based nonlinear solver that uses a line search
17132b4ed282SShri Abhyankar 
17142b4ed282SShri Abhyankar    Options Database:
1715be6adb11SBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1716be6adb11SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1717be6adb11SBarry Smith .   -snes_ls_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)
1718be6adb11SBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1719be6adb11SBarry Smith -   -snes_ls_monitor - print information about progress of line searches
17202b4ed282SShri Abhyankar 
17212b4ed282SShri Abhyankar 
17222b4ed282SShri Abhyankar    Level: beginner
17232b4ed282SShri Abhyankar 
17242b4ed282SShri Abhyankar .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
17252b4ed282SShri Abhyankar            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
17262b4ed282SShri Abhyankar            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
17272b4ed282SShri Abhyankar 
17282b4ed282SShri Abhyankar M*/
17292b4ed282SShri Abhyankar EXTERN_C_BEGIN
17302b4ed282SShri Abhyankar #undef __FUNCT__
17312b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI"
17327087cfbeSBarry Smith PetscErrorCode  SNESCreate_VI(SNES snes)
17332b4ed282SShri Abhyankar {
17342b4ed282SShri Abhyankar   PetscErrorCode ierr;
17352b4ed282SShri Abhyankar   SNES_VI      *vi;
17362b4ed282SShri Abhyankar 
17372b4ed282SShri Abhyankar   PetscFunctionBegin;
17382b4ed282SShri Abhyankar   snes->ops->setup	     = SNESSetUp_VI;
173969c03620SShri Abhyankar   snes->ops->solve	     = SNESSolveVI_SS;
17402b4ed282SShri Abhyankar   snes->ops->destroy	     = SNESDestroy_VI;
17412b4ed282SShri Abhyankar   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
17422b4ed282SShri Abhyankar   snes->ops->view            = SNESView_VI;
17432b4ed282SShri Abhyankar   snes->ops->converged       = SNESDefaultConverged_VI;
17442b4ed282SShri Abhyankar 
17452b4ed282SShri Abhyankar   ierr                   = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
17462b4ed282SShri Abhyankar   snes->data    	 = (void*)vi;
17472b4ed282SShri Abhyankar   vi->alpha		 = 1.e-4;
17482b4ed282SShri Abhyankar   vi->maxstep		 = 1.e8;
17492b4ed282SShri Abhyankar   vi->minlambda         = 1.e-12;
17507fe79bc4SShri Abhyankar   vi->LineSearch        = SNESLineSearchNo_VI;
17512b4ed282SShri Abhyankar   vi->lsP               = PETSC_NULL;
17522b4ed282SShri Abhyankar   vi->postcheckstep     = PETSC_NULL;
17532b4ed282SShri Abhyankar   vi->postcheck         = PETSC_NULL;
17542b4ed282SShri Abhyankar   vi->precheckstep      = PETSC_NULL;
17552b4ed282SShri Abhyankar   vi->precheck          = PETSC_NULL;
17562b4ed282SShri Abhyankar   vi->const_tol         =  2.2204460492503131e-16;
17572b4ed282SShri Abhyankar 
17582b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
17592b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
17602b4ed282SShri Abhyankar 
17612b4ed282SShri Abhyankar   PetscFunctionReturn(0);
17622b4ed282SShri Abhyankar }
17632b4ed282SShri Abhyankar EXTERN_C_END
1764