xref: /petsc/src/snes/impls/vi/ss/viss.c (revision fbda97447efb739befdf7134ad897ba9ea04a9a4)
1c2fc9fa9SBarry Smith 
2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/
3c2fc9fa9SBarry Smith 
4f6dfbefdSBarry Smith /*@
5c2fc9fa9SBarry Smith   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
6c2fc9fa9SBarry Smith 
7c2fc9fa9SBarry Smith    Input Parameter:
8f6dfbefdSBarry Smith .  phi - the `Vec` holding the evaluation of the semismooth function
9c2fc9fa9SBarry Smith 
10f6dfbefdSBarry Smith    Output Parameters:
11f6dfbefdSBarry Smith +  merit - the merit function 1/2 ||phi||^2
12f6dfbefdSBarry Smith -  phinorm -  the two-norm of the vector, ||phi||
13c2fc9fa9SBarry Smith 
14*fbda9744SBarry Smith    Level: developer
15*fbda9744SBarry Smith 
16f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeFunction()`
17f6dfbefdSBarry Smith @*/
18f6dfbefdSBarry Smith PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm) {
19c2fc9fa9SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(phi, NORM_2, phinorm));
219566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(phi, NORM_2, phinorm));
22c2fc9fa9SBarry Smith   *merit = 0.5 * (*phinorm) * (*phinorm);
23c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
24c2fc9fa9SBarry Smith }
25c2fc9fa9SBarry Smith 
269371c9d4SSatish Balay static inline PetscScalar Phi(PetscScalar a, PetscScalar b) {
27c2fc9fa9SBarry Smith   return a + b - PetscSqrtScalar(a * a + b * b);
28c2fc9fa9SBarry Smith }
29c2fc9fa9SBarry Smith 
309371c9d4SSatish Balay static inline PetscScalar DPhi(PetscScalar a, PetscScalar b) {
31c2fc9fa9SBarry Smith   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b);
32c2fc9fa9SBarry Smith   else return .5;
33c2fc9fa9SBarry Smith }
34c2fc9fa9SBarry Smith 
35f6dfbefdSBarry Smith /*@
36f6dfbefdSBarry Smith    SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear
37f6dfbefdSBarry Smith    equations in semismooth form.
38c2fc9fa9SBarry Smith 
39c2fc9fa9SBarry Smith    Input Parameters:
40f6dfbefdSBarry Smith +  snes - the SNES context
41d5f1b7e6SEd Bueler .  X - current iterate
42f6dfbefdSBarry Smith -  functx - user defined function context
43c2fc9fa9SBarry Smith 
44f6dfbefdSBarry Smith    Output Parameter:
45f6dfbefdSBarry Smith .  phi - the evaluation of Semismooth function at X
46c2fc9fa9SBarry Smith 
47*fbda9744SBarry Smith    Level: developer
48*fbda9744SBarry Smith 
49f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()`
50f6dfbefdSBarry Smith @*/
51f6dfbefdSBarry Smith PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx) {
52f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
53c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu, F = snes->vec_func;
545edff71fSBarry Smith   PetscScalar       *phi_arr, *f_arr, *l, *u;
555edff71fSBarry Smith   const PetscScalar *x_arr;
56c2fc9fa9SBarry Smith   PetscInt           i, nlocal;
57c2fc9fa9SBarry Smith 
58c2fc9fa9SBarry Smith   PetscFunctionBegin;
599566063dSJacob Faibussowitsch   PetscCall((*vi->computeuserfunction)(snes, X, F, functx));
609566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x_arr));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f_arr));
639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xl, &l));
649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xu, &u));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(phi, &phi_arr));
66c2fc9fa9SBarry Smith 
67c2fc9fa9SBarry Smith   for (i = 0; i < nlocal; i++) {
68e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
69c2fc9fa9SBarry Smith       phi_arr[i] = f_arr[i];
70e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
71c2fc9fa9SBarry Smith       phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]);
72e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
73c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]);
74c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) {
75c2fc9fa9SBarry Smith       phi_arr[i] = l[i] - x_arr[i];
76c2fc9fa9SBarry Smith     } else { /* both bounds on variable */
77c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i]));
78c2fc9fa9SBarry Smith     }
79c2fc9fa9SBarry Smith   }
80c2fc9fa9SBarry Smith 
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x_arr));
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f_arr));
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xl, &l));
849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xu, &u));
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(phi, &phi_arr));
86c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
87c2fc9fa9SBarry Smith }
88c2fc9fa9SBarry Smith 
89c2fc9fa9SBarry Smith /*
90c2fc9fa9SBarry Smith    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
91c2fc9fa9SBarry Smith                                           the semismooth jacobian.
92c2fc9fa9SBarry Smith */
939371c9d4SSatish Balay PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db) {
94c2fc9fa9SBarry Smith   PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2;
95c2fc9fa9SBarry Smith   PetscInt     i, nlocal;
96c2fc9fa9SBarry Smith 
97c2fc9fa9SBarry Smith   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
1009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(snes->xl, &l));
1019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(snes->xu, &u));
1029566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Da, &da));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Db, &db));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
105c2fc9fa9SBarry Smith 
106c2fc9fa9SBarry Smith   for (i = 0; i < nlocal; i++) {
107e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
108c2fc9fa9SBarry Smith       da[i] = 0;
109c2fc9fa9SBarry Smith       db[i] = 1;
110e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
111c2fc9fa9SBarry Smith       da[i] = DPhi(u[i] - x[i], -f[i]);
112c2fc9fa9SBarry Smith       db[i] = DPhi(-f[i], u[i] - x[i]);
113e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
114c2fc9fa9SBarry Smith       da[i] = DPhi(x[i] - l[i], f[i]);
115c2fc9fa9SBarry Smith       db[i] = DPhi(f[i], x[i] - l[i]);
116c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) { /* fixed variable */
117c2fc9fa9SBarry Smith       da[i] = 1;
118c2fc9fa9SBarry Smith       db[i] = 0;
119c2fc9fa9SBarry Smith     } else { /* upper and lower bounds on variable */
120c2fc9fa9SBarry Smith       da1   = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
121c2fc9fa9SBarry Smith       db1   = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]);
122c2fc9fa9SBarry Smith       da2   = DPhi(u[i] - x[i], -f[i]);
123c2fc9fa9SBarry Smith       db2   = DPhi(-f[i], u[i] - x[i]);
124c2fc9fa9SBarry Smith       da[i] = da1 + db1 * da2;
125c2fc9fa9SBarry Smith       db[i] = db1 * db2;
126c2fc9fa9SBarry Smith     }
127c2fc9fa9SBarry Smith   }
128c2fc9fa9SBarry Smith 
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(snes->xl, &l));
1329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(snes->xu, &u));
1339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Da, &da));
1349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Db, &db));
135c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
136c2fc9fa9SBarry Smith }
137c2fc9fa9SBarry Smith 
138c2fc9fa9SBarry Smith /*
139c2fc9fa9SBarry Smith    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.
140c2fc9fa9SBarry Smith 
141c2fc9fa9SBarry Smith    Input Parameters:
142c2fc9fa9SBarry Smith .  Da       - Diagonal shift vector for the semismooth jacobian.
143c2fc9fa9SBarry Smith .  Db       - Row scaling vector for the semismooth jacobian.
144c2fc9fa9SBarry Smith 
145c2fc9fa9SBarry Smith    Output Parameters:
146c2fc9fa9SBarry Smith .  jac      - semismooth jacobian
147c2fc9fa9SBarry Smith .  jac_pre  - optional preconditioning matrix
148c2fc9fa9SBarry Smith 
149f6dfbefdSBarry Smith    Note:
150c2fc9fa9SBarry Smith    The semismooth jacobian matrix is given by
151c2fc9fa9SBarry Smith    jac = Da + Db*jacfun
152c2fc9fa9SBarry Smith    where Db is the row scaling matrix stored as a vector,
153c2fc9fa9SBarry Smith          Da is the diagonal perturbation matrix stored as a vector
154c2fc9fa9SBarry Smith    and   jacfun is the jacobian of the original nonlinear function.
155c2fc9fa9SBarry Smith */
1569371c9d4SSatish Balay PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db) {
157c2fc9fa9SBarry Smith   /* Do row scaling  and add diagonal perturbation */
158362febeeSStefano Zampini   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(jac, Db, NULL));
1609566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES));
161c2fc9fa9SBarry Smith   if (jac != jac_pre) { /* If jac and jac_pre are different */
1629566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(jac_pre, Db, NULL));
1639566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES));
164c2fc9fa9SBarry Smith   }
165c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
166c2fc9fa9SBarry Smith }
167c2fc9fa9SBarry Smith 
168c2fc9fa9SBarry Smith /*
169c2fc9fa9SBarry Smith    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
170c2fc9fa9SBarry Smith 
171c2fc9fa9SBarry Smith    Input Parameters:
172c2fc9fa9SBarry Smith    phi - semismooth function.
173c2fc9fa9SBarry Smith    H   - semismooth jacobian
174c2fc9fa9SBarry Smith 
175f6dfbefdSBarry Smith    Output Parameter:
176c2fc9fa9SBarry Smith    dpsi - merit function gradient
177c2fc9fa9SBarry Smith 
178f6dfbefdSBarry Smith    Note:
179c2fc9fa9SBarry Smith   The merit function gradient is computed as follows
180c2fc9fa9SBarry Smith         dpsi = H^T*phi
181c2fc9fa9SBarry Smith */
1829371c9d4SSatish Balay PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) {
183c2fc9fa9SBarry Smith   PetscFunctionBegin;
1849566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(H, phi, dpsi));
185c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
186c2fc9fa9SBarry Smith }
187c2fc9fa9SBarry Smith 
188c2fc9fa9SBarry Smith /*
189f450aa47SBarry Smith    SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton
190c2fc9fa9SBarry Smith    method using a line search.
191c2fc9fa9SBarry Smith 
192f6dfbefdSBarry Smith    Input Parameter:
193c2fc9fa9SBarry Smith .  snes - the SNES context
194c2fc9fa9SBarry Smith 
195c2fc9fa9SBarry Smith    Application Interface Routine: SNESSolve()
196c2fc9fa9SBarry Smith 
197f6dfbefdSBarry Smith    Note:
198c2fc9fa9SBarry Smith    This implements essentially a semismooth Newton method with a
199d5f1b7e6SEd Bueler    line search. The default line search does not do any line search
200d5f1b7e6SEd Bueler    but rather takes a full Newton step.
201016d19f9SBarry Smith 
20204d7464bSBarry Smith    Developer Note: the code in this file should be slightly modified so that this routine need not exist and the SNESSolve_NEWTONLS() routine is called directly with the appropriate wrapped function and Jacobian evaluations
203016d19f9SBarry Smith 
204c2fc9fa9SBarry Smith */
2059371c9d4SSatish Balay PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) {
206f450aa47SBarry Smith   SNES_VINEWTONSSLS   *vi = (SNES_VINEWTONSSLS *)snes->data;
207c2fc9fa9SBarry Smith   PetscInt             maxits, i, lits;
208422a814eSBarry Smith   SNESLineSearchReason lssucceed;
209c2fc9fa9SBarry Smith   PetscReal            gnorm, xnorm = 0, ynorm;
2109bd66eb0SPeter Brune   Vec                  Y, X, F;
211c2fc9fa9SBarry Smith   KSPConvergedReason   kspreason;
2126cab3a1bSJed Brown   DM                   dm;
213942e3340SBarry Smith   DMSNES               sdm;
214c2fc9fa9SBarry Smith 
215c2fc9fa9SBarry Smith   PetscFunctionBegin;
2169566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
2179566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
2181aa26658SKarl Rupp 
21922c6f798SBarry Smith   vi->computeuserfunction   = sdm->ops->computefunction;
22022c6f798SBarry Smith   sdm->ops->computefunction = SNESVIComputeFunction;
221c2fc9fa9SBarry Smith 
222c2fc9fa9SBarry Smith   snes->numFailures            = 0;
223c2fc9fa9SBarry Smith   snes->numLinearSolveFailures = 0;
224c2fc9fa9SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
225c2fc9fa9SBarry Smith 
226c2fc9fa9SBarry Smith   maxits = snes->max_its;  /* maximum number of iterations */
227c2fc9fa9SBarry Smith   X      = snes->vec_sol;  /* solution vector */
228c2fc9fa9SBarry Smith   F      = snes->vec_func; /* residual vector */
229c2fc9fa9SBarry Smith   Y      = snes->work[0];  /* work vectors */
230c2fc9fa9SBarry Smith 
2319566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
232c2fc9fa9SBarry Smith   snes->iter = 0;
233c2fc9fa9SBarry Smith   snes->norm = 0.0;
2349566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
235c2fc9fa9SBarry Smith 
2369566063dSJacob Faibussowitsch   PetscCall(SNESVIProjectOntoBounds(snes, X));
2379566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, vi->phi));
238c2fc9fa9SBarry Smith   if (snes->domainerror) {
239c2fc9fa9SBarry Smith     snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
24022c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
241c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
242c2fc9fa9SBarry Smith   }
243c2fc9fa9SBarry Smith   /* Compute Merit function */
2449566063dSJacob Faibussowitsch   PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm));
245c2fc9fa9SBarry Smith 
2469566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x||  */
2479566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
248422a814eSBarry Smith   SNESCheckFunctionNorm(snes, vi->merit);
249c2fc9fa9SBarry Smith 
2509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
251c2fc9fa9SBarry Smith   snes->norm = vi->phinorm;
2529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2539566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0));
2549566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, vi->phinorm));
255c2fc9fa9SBarry Smith 
256c2fc9fa9SBarry Smith   /* test convergence */
257dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, vi->phinorm, &snes->reason, snes->cnvP);
258c2fc9fa9SBarry Smith   if (snes->reason) {
25922c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
260c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
261c2fc9fa9SBarry Smith   }
262c2fc9fa9SBarry Smith 
263c2fc9fa9SBarry Smith   for (i = 0; i < maxits; i++) {
264c2fc9fa9SBarry Smith     /* Call general purpose update function */
265dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
266c2fc9fa9SBarry Smith 
267c2fc9fa9SBarry Smith     /* Solve J Y = Phi, where J is the semismooth jacobian */
268b9600128SPeter Brune 
269b9600128SPeter Brune     /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/
27022c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
2719566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
27207b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
27322c6f798SBarry Smith     sdm->ops->computefunction = SNESVIComputeFunction;
274b9600128SPeter Brune 
275c2fc9fa9SBarry Smith     /* Get the diagonal shift and row scaling vectors */
2769566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db));
277c2fc9fa9SBarry Smith     /* Compute the semismooth jacobian */
2789566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db));
279c2fc9fa9SBarry Smith     /* Compute the merit function gradient */
2809566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi));
2819566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2829566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, vi->phi, Y));
2839566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason));
284c2fc9fa9SBarry Smith 
285c2fc9fa9SBarry Smith     if (kspreason < 0) {
286c2fc9fa9SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
28763a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures));
288c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
289c2fc9fa9SBarry Smith         break;
290c2fc9fa9SBarry Smith       }
291c2fc9fa9SBarry Smith     }
2929566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
293c2fc9fa9SBarry Smith     snes->linear_its += lits;
29463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
295c2fc9fa9SBarry Smith     /*
2966b2b7091SBarry Smith     if (snes->ops->precheck) {
297c2fc9fa9SBarry Smith       PetscBool changed_y = PETSC_FALSE;
298dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y);
299c2fc9fa9SBarry Smith     }
300c2fc9fa9SBarry Smith 
3011baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W));
302c2fc9fa9SBarry Smith     */
303c2fc9fa9SBarry Smith     /* Compute a (scaled) negative update in the line search routine:
304c2fc9fa9SBarry Smith          Y <- X - lambda*Y
305c2fc9fa9SBarry Smith        and evaluate G = function(Y) (depends on the line search).
306c2fc9fa9SBarry Smith     */
3079566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, snes->vec_sol_update));
3089371c9d4SSatish Balay     ynorm = 1;
3099371c9d4SSatish Balay     gnorm = vi->phinorm;
3109566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y));
3119566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
3129566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm));
3139566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)vi->phinorm, (double)gnorm, (double)ynorm, (int)lssucceed));
314c2fc9fa9SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
315c2fc9fa9SBarry Smith     if (snes->domainerror) {
316c2fc9fa9SBarry Smith       snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
31722c6f798SBarry Smith       sdm->ops->computefunction = vi->computeuserfunction;
318c2fc9fa9SBarry Smith       PetscFunctionReturn(0);
319c2fc9fa9SBarry Smith     }
320422a814eSBarry Smith     if (lssucceed) {
321c2fc9fa9SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
322c2fc9fa9SBarry Smith         PetscBool ismin;
323c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3249566063dSJacob Faibussowitsch         PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin));
325c2fc9fa9SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
326c2fc9fa9SBarry Smith         break;
327c2fc9fa9SBarry Smith       }
328c2fc9fa9SBarry Smith     }
329c2fc9fa9SBarry Smith     /* Update function and solution vectors */
330c2fc9fa9SBarry Smith     vi->phinorm = gnorm;
331c2fc9fa9SBarry Smith     vi->merit   = 0.5 * vi->phinorm * vi->phinorm;
332c2fc9fa9SBarry Smith     /* Monitor convergence */
3339566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
334c2fc9fa9SBarry Smith     snes->iter  = i + 1;
335c2fc9fa9SBarry Smith     snes->norm  = vi->phinorm;
336c1e67a49SFande Kong     snes->xnorm = xnorm;
337c1e67a49SFande Kong     snes->ynorm = ynorm;
3389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3399566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
3409566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
341c2fc9fa9SBarry Smith     /* Test for convergence, xnorm = || X || */
3429566063dSJacob Faibussowitsch     if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
343dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, vi->phinorm, &snes->reason, snes->cnvP);
344c2fc9fa9SBarry Smith     if (snes->reason) break;
345c2fc9fa9SBarry Smith   }
346c2fc9fa9SBarry Smith   if (i == maxits) {
34763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
348c2fc9fa9SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
349c2fc9fa9SBarry Smith   }
35022c6f798SBarry Smith   sdm->ops->computefunction = vi->computeuserfunction;
351c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
352c2fc9fa9SBarry Smith }
353c2fc9fa9SBarry Smith 
354c2fc9fa9SBarry Smith /*
355f450aa47SBarry Smith    SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use
356c2fc9fa9SBarry Smith    of the SNES nonlinear solver.
357c2fc9fa9SBarry Smith 
358c2fc9fa9SBarry Smith    Input Parameter:
359c2fc9fa9SBarry Smith .  snes - the SNES context
360c2fc9fa9SBarry Smith 
361c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetUp()
362c2fc9fa9SBarry Smith 
363f6dfbefdSBarry Smith    Note:
364c2fc9fa9SBarry Smith    For basic use of the SNES solvers, the user need not explicitly call
365c2fc9fa9SBarry Smith    SNESSetUp(), since these actions will automatically occur during
366c2fc9fa9SBarry Smith    the call to SNESSolve().
367c2fc9fa9SBarry Smith  */
3689371c9d4SSatish Balay PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) {
369f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
370c2fc9fa9SBarry Smith 
371c2fc9fa9SBarry Smith   PetscFunctionBegin;
3729566063dSJacob Faibussowitsch   PetscCall(SNESSetUp_VI(snes));
3739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi));
3749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->phi));
3759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->Da));
3769566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->Db));
3779566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->z));
3789566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->t));
379c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
380c2fc9fa9SBarry Smith }
381f6dfbefdSBarry Smith 
3829371c9d4SSatish Balay PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) {
383f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
384c2fc9fa9SBarry Smith 
385c2fc9fa9SBarry Smith   PetscFunctionBegin;
3869566063dSJacob Faibussowitsch   PetscCall(SNESReset_VI(snes));
3879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->dpsi));
3889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->phi));
3899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->Da));
3909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->Db));
3919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->z));
3929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->t));
393c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
394c2fc9fa9SBarry Smith }
395c2fc9fa9SBarry Smith 
396c2fc9fa9SBarry Smith /*
397f450aa47SBarry Smith    SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method.
398c2fc9fa9SBarry Smith 
399c2fc9fa9SBarry Smith    Input Parameter:
400c2fc9fa9SBarry Smith .  snes - the SNES context
401c2fc9fa9SBarry Smith 
402c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetFromOptions()
403c2fc9fa9SBarry Smith */
4049371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject) {
405c2fc9fa9SBarry Smith   PetscFunctionBegin;
406dbbe0bcdSBarry Smith   PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject));
407d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options");
408d0609cedSBarry Smith   PetscOptionsHeadEnd();
409c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
410c2fc9fa9SBarry Smith }
411c2fc9fa9SBarry Smith 
412c2fc9fa9SBarry Smith /*MC
413f450aa47SBarry Smith       SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method
414c2fc9fa9SBarry Smith 
415f6dfbefdSBarry Smith    Options Database Keys:
416d5f1b7e6SEd Bueler +   -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method
41761589011SJed Brown -   -snes_vi_monitor - prints the number of active constraints at each iteration.
41861589011SJed Brown 
419c2fc9fa9SBarry Smith    Level: beginner
420c2fc9fa9SBarry Smith 
421b80f3ac1SShri Abhyankar    References:
422606c0280SSatish Balay +  * -  T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth
423b80f3ac1SShri Abhyankar      algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001).
424606c0280SSatish Balay -  * -  T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale
425d5f1b7e6SEd Bueler      Applications, Optimization Methods and Software, 21 (2006).
426b80f3ac1SShri Abhyankar 
427f6dfbefdSBarry Smith    Notes:
428f6dfbefdSBarry Smith    This family of algorithm is much like an interior point method.
429c2fc9fa9SBarry Smith 
430f6dfbefdSBarry Smith    The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems
431f6dfbefdSBarry Smith 
432f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
433c2fc9fa9SBarry Smith M*/
4349371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) {
435f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi;
436d8d34be6SBarry Smith   SNESLineSearch     linesearch;
437c2fc9fa9SBarry Smith 
438c2fc9fa9SBarry Smith   PetscFunctionBegin;
439f450aa47SBarry Smith   snes->ops->reset          = SNESReset_VINEWTONSSLS;
440f450aa47SBarry Smith   snes->ops->setup          = SNESSetUp_VINEWTONSSLS;
441f450aa47SBarry Smith   snes->ops->solve          = SNESSolve_VINEWTONSSLS;
442c2fc9fa9SBarry Smith   snes->ops->destroy        = SNESDestroy_VI;
443f450aa47SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS;
4440298fd71SBarry Smith   snes->ops->view           = NULL;
445c2fc9fa9SBarry Smith 
446c2fc9fa9SBarry Smith   snes->usesksp = PETSC_TRUE;
447efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
448c2fc9fa9SBarry Smith 
4499566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
450ec786807SJed Brown   if (!((PetscObject)linesearch)->type_name) {
4519566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
4529566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0));
453ec786807SJed Brown   }
454d8d34be6SBarry Smith 
4554fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
4564fc747eaSLawrence Mitchell 
4579566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &vi));
458c2fc9fa9SBarry Smith   snes->data = (void *)vi;
459c2fc9fa9SBarry Smith 
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI));
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI));
462c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
463c2fc9fa9SBarry Smith }
464