xref: /petsc/src/snes/impls/vi/ss/viss.c (revision ceaaa4989964adb3f5eb130cb04b8f49c83e49c9)
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 
14fbda9744SBarry Smith   Level: developer
15fbda9744SBarry Smith 
16f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeFunction()`
17f6dfbefdSBarry Smith @*/
18d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm)
19d71ae5a4SJacob Faibussowitsch {
20c2fc9fa9SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(phi, NORM_2, phinorm));
229566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(phi, NORM_2, phinorm));
23c2fc9fa9SBarry Smith   *merit = 0.5 * (*phinorm) * (*phinorm);
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25c2fc9fa9SBarry Smith }
26c2fc9fa9SBarry Smith 
27d71ae5a4SJacob Faibussowitsch static inline PetscScalar Phi(PetscScalar a, PetscScalar b)
28d71ae5a4SJacob Faibussowitsch {
29c2fc9fa9SBarry Smith   return a + b - PetscSqrtScalar(a * a + b * b);
30c2fc9fa9SBarry Smith }
31c2fc9fa9SBarry Smith 
32d71ae5a4SJacob Faibussowitsch static inline PetscScalar DPhi(PetscScalar a, PetscScalar b)
33d71ae5a4SJacob Faibussowitsch {
34c2fc9fa9SBarry Smith   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b);
35c2fc9fa9SBarry Smith   else return .5;
36c2fc9fa9SBarry Smith }
37c2fc9fa9SBarry Smith 
38f6dfbefdSBarry Smith /*@
39f6dfbefdSBarry Smith   SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear
40f6dfbefdSBarry Smith   equations in semismooth form.
41c2fc9fa9SBarry Smith 
42c2fc9fa9SBarry Smith   Input Parameters:
43f6dfbefdSBarry Smith + snes   - the SNES context
44d5f1b7e6SEd Bueler . X      - current iterate
45f6dfbefdSBarry Smith - functx - user defined function context
46c2fc9fa9SBarry Smith 
47f6dfbefdSBarry Smith   Output Parameter:
48f6dfbefdSBarry Smith . phi - the evaluation of Semismooth function at X
49c2fc9fa9SBarry Smith 
50fbda9744SBarry Smith   Level: developer
51fbda9744SBarry Smith 
52f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()`
53f6dfbefdSBarry Smith @*/
54d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx)
55d71ae5a4SJacob Faibussowitsch {
56f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
57c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu, F = snes->vec_func;
585edff71fSBarry Smith   PetscScalar       *phi_arr, *f_arr, *l, *u;
595edff71fSBarry Smith   const PetscScalar *x_arr;
60c2fc9fa9SBarry Smith   PetscInt           i, nlocal;
61c2fc9fa9SBarry Smith 
62c2fc9fa9SBarry Smith   PetscFunctionBegin;
639566063dSJacob Faibussowitsch   PetscCall((*vi->computeuserfunction)(snes, X, F, functx));
649566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x_arr));
669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f_arr));
679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xl, &l));
689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xu, &u));
699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(phi, &phi_arr));
70c2fc9fa9SBarry Smith 
71c2fc9fa9SBarry Smith   for (i = 0; i < nlocal; i++) {
72e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
73c2fc9fa9SBarry Smith       phi_arr[i] = f_arr[i];
74e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
75c2fc9fa9SBarry Smith       phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]);
76e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
77c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]);
78c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) {
79c2fc9fa9SBarry Smith       phi_arr[i] = l[i] - x_arr[i];
80c2fc9fa9SBarry Smith     } else { /* both bounds on variable */
81c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i]));
82c2fc9fa9SBarry Smith     }
83c2fc9fa9SBarry Smith   }
84c2fc9fa9SBarry Smith 
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x_arr));
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f_arr));
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xl, &l));
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xu, &u));
899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(phi, &phi_arr));
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91c2fc9fa9SBarry Smith }
92c2fc9fa9SBarry Smith 
93c2fc9fa9SBarry Smith /*
94c2fc9fa9SBarry Smith    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
95c2fc9fa9SBarry Smith                                           the semismooth jacobian.
96c2fc9fa9SBarry Smith */
97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db)
98d71ae5a4SJacob Faibussowitsch {
99c2fc9fa9SBarry Smith   PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2;
100c2fc9fa9SBarry Smith   PetscInt     i, nlocal;
101c2fc9fa9SBarry Smith 
102c2fc9fa9SBarry Smith   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(snes->xl, &l));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(snes->xu, &u));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Da, &da));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Db, &db));
1099566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
110c2fc9fa9SBarry Smith 
111c2fc9fa9SBarry Smith   for (i = 0; i < nlocal; i++) {
112e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
113c2fc9fa9SBarry Smith       da[i] = 0;
114c2fc9fa9SBarry Smith       db[i] = 1;
115e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
116c2fc9fa9SBarry Smith       da[i] = DPhi(u[i] - x[i], -f[i]);
117c2fc9fa9SBarry Smith       db[i] = DPhi(-f[i], u[i] - x[i]);
118e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
119c2fc9fa9SBarry Smith       da[i] = DPhi(x[i] - l[i], f[i]);
120c2fc9fa9SBarry Smith       db[i] = DPhi(f[i], x[i] - l[i]);
121c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) { /* fixed variable */
122c2fc9fa9SBarry Smith       da[i] = 1;
123c2fc9fa9SBarry Smith       db[i] = 0;
124c2fc9fa9SBarry Smith     } else { /* upper and lower bounds on variable */
125c2fc9fa9SBarry Smith       da1   = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
126c2fc9fa9SBarry Smith       db1   = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]);
127c2fc9fa9SBarry Smith       da2   = DPhi(u[i] - x[i], -f[i]);
128c2fc9fa9SBarry Smith       db2   = DPhi(-f[i], u[i] - x[i]);
129c2fc9fa9SBarry Smith       da[i] = da1 + db1 * da2;
130c2fc9fa9SBarry Smith       db[i] = db1 * db2;
131c2fc9fa9SBarry Smith     }
132c2fc9fa9SBarry Smith   }
133c2fc9fa9SBarry Smith 
1349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(snes->xl, &l));
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(snes->xu, &u));
1389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Da, &da));
1399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Db, &db));
1403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141c2fc9fa9SBarry Smith }
142c2fc9fa9SBarry Smith 
143c2fc9fa9SBarry Smith /*
144c2fc9fa9SBarry 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.
145c2fc9fa9SBarry Smith 
146c2fc9fa9SBarry Smith    Input Parameters:
147c2fc9fa9SBarry Smith .  Da       - Diagonal shift vector for the semismooth jacobian.
148c2fc9fa9SBarry Smith .  Db       - Row scaling vector for the semismooth jacobian.
149c2fc9fa9SBarry Smith 
150c2fc9fa9SBarry Smith    Output Parameters:
151c2fc9fa9SBarry Smith .  jac      - semismooth jacobian
152c2fc9fa9SBarry Smith .  jac_pre  - optional preconditioning matrix
153c2fc9fa9SBarry Smith 
154f6dfbefdSBarry Smith    Note:
155c2fc9fa9SBarry Smith    The semismooth jacobian matrix is given by
156c2fc9fa9SBarry Smith    jac = Da + Db*jacfun
157c2fc9fa9SBarry Smith    where Db is the row scaling matrix stored as a vector,
158c2fc9fa9SBarry Smith          Da is the diagonal perturbation matrix stored as a vector
159c2fc9fa9SBarry Smith    and   jacfun is the jacobian of the original nonlinear function.
160c2fc9fa9SBarry Smith */
161d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db)
162d71ae5a4SJacob Faibussowitsch {
163c2fc9fa9SBarry Smith   /* Do row scaling  and add diagonal perturbation */
164362febeeSStefano Zampini   PetscFunctionBegin;
1659566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(jac, Db, NULL));
1669566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES));
167c2fc9fa9SBarry Smith   if (jac != jac_pre) { /* If jac and jac_pre are different */
1689566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(jac_pre, Db, NULL));
1699566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES));
170c2fc9fa9SBarry Smith   }
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172c2fc9fa9SBarry Smith }
173c2fc9fa9SBarry Smith 
174*ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars
175c2fc9fa9SBarry Smith /*
176c2fc9fa9SBarry Smith   SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
177c2fc9fa9SBarry Smith 
178c2fc9fa9SBarry Smith   Input Parameters:
179c2fc9fa9SBarry Smith    phi - semismooth function.
180c2fc9fa9SBarry Smith    H   - semismooth jacobian
181c2fc9fa9SBarry Smith 
182f6dfbefdSBarry Smith   Output Parameter:
183c2fc9fa9SBarry Smith    dpsi - merit function gradient
184c2fc9fa9SBarry Smith 
185f6dfbefdSBarry Smith   Note:
186c2fc9fa9SBarry Smith   The merit function gradient is computed as follows
187c2fc9fa9SBarry Smith   dpsi = H^T*phi
188c2fc9fa9SBarry Smith */
189*ceaaa498SBarry Smith static PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
190d71ae5a4SJacob Faibussowitsch {
191c2fc9fa9SBarry Smith   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(H, phi, dpsi));
1933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194c2fc9fa9SBarry Smith }
195c2fc9fa9SBarry Smith 
196d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes)
197d71ae5a4SJacob Faibussowitsch {
198f450aa47SBarry Smith   SNES_VINEWTONSSLS   *vi = (SNES_VINEWTONSSLS *)snes->data;
199c2fc9fa9SBarry Smith   PetscInt             maxits, i, lits;
200422a814eSBarry Smith   SNESLineSearchReason lssucceed;
201c2fc9fa9SBarry Smith   PetscReal            gnorm, xnorm = 0, ynorm;
2029bd66eb0SPeter Brune   Vec                  Y, X, F;
203c2fc9fa9SBarry Smith   KSPConvergedReason   kspreason;
2046cab3a1bSJed Brown   DM                   dm;
205942e3340SBarry Smith   DMSNES               sdm;
206c2fc9fa9SBarry Smith 
207c2fc9fa9SBarry Smith   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
2099566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
2101aa26658SKarl Rupp 
21122c6f798SBarry Smith   vi->computeuserfunction   = sdm->ops->computefunction;
21222c6f798SBarry Smith   sdm->ops->computefunction = SNESVIComputeFunction;
213c2fc9fa9SBarry Smith 
214c2fc9fa9SBarry Smith   snes->numFailures            = 0;
215c2fc9fa9SBarry Smith   snes->numLinearSolveFailures = 0;
216c2fc9fa9SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
217c2fc9fa9SBarry Smith 
218c2fc9fa9SBarry Smith   maxits = snes->max_its;  /* maximum number of iterations */
219c2fc9fa9SBarry Smith   X      = snes->vec_sol;  /* solution vector */
220c2fc9fa9SBarry Smith   F      = snes->vec_func; /* residual vector */
221c2fc9fa9SBarry Smith   Y      = snes->work[0];  /* work vectors */
222c2fc9fa9SBarry Smith 
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
224c2fc9fa9SBarry Smith   snes->iter = 0;
225c2fc9fa9SBarry Smith   snes->norm = 0.0;
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
227c2fc9fa9SBarry Smith 
2289566063dSJacob Faibussowitsch   PetscCall(SNESVIProjectOntoBounds(snes, X));
2299566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, vi->phi));
230c2fc9fa9SBarry Smith   if (snes->domainerror) {
231c2fc9fa9SBarry Smith     snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
23222c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
2333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
234c2fc9fa9SBarry Smith   }
235c2fc9fa9SBarry Smith   /* Compute Merit function */
2369566063dSJacob Faibussowitsch   PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm));
237c2fc9fa9SBarry Smith 
2389566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x||  */
2399566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
240422a814eSBarry Smith   SNESCheckFunctionNorm(snes, vi->merit);
241c2fc9fa9SBarry Smith 
2429566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
243c2fc9fa9SBarry Smith   snes->norm = vi->phinorm;
2449566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2459566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0));
246c2fc9fa9SBarry Smith 
247c2fc9fa9SBarry Smith   /* test convergence */
2482d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, vi->phinorm));
2492d157150SStefano Zampini   PetscCall(SNESMonitor(snes, 0, vi->phinorm));
250c2fc9fa9SBarry Smith   if (snes->reason) {
25122c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
2523ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
253c2fc9fa9SBarry Smith   }
254c2fc9fa9SBarry Smith 
255c2fc9fa9SBarry Smith   for (i = 0; i < maxits; i++) {
256c2fc9fa9SBarry Smith     /* Call general purpose update function */
257dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
258c2fc9fa9SBarry Smith 
259c2fc9fa9SBarry Smith     /* Solve J Y = Phi, where J is the semismooth jacobian */
260b9600128SPeter Brune 
261b9600128SPeter Brune     /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/
26222c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
2639566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
26407b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
26522c6f798SBarry Smith     sdm->ops->computefunction = SNESVIComputeFunction;
266b9600128SPeter Brune 
267c2fc9fa9SBarry Smith     /* Get the diagonal shift and row scaling vectors */
2689566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db));
269c2fc9fa9SBarry Smith     /* Compute the semismooth jacobian */
2709566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db));
271c2fc9fa9SBarry Smith     /* Compute the merit function gradient */
2729566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi));
2739566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2749566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, vi->phi, Y));
2759566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason));
276c2fc9fa9SBarry Smith 
277c2fc9fa9SBarry Smith     if (kspreason < 0) {
278c2fc9fa9SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
27963a3b9bcSJacob 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));
280c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
281c2fc9fa9SBarry Smith         break;
282c2fc9fa9SBarry Smith       }
283c2fc9fa9SBarry Smith     }
2849566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
285c2fc9fa9SBarry Smith     snes->linear_its += lits;
28663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
287c2fc9fa9SBarry Smith     /*
2886b2b7091SBarry Smith     if (snes->ops->precheck) {
289c2fc9fa9SBarry Smith       PetscBool changed_y = PETSC_FALSE;
290dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y);
291c2fc9fa9SBarry Smith     }
292c2fc9fa9SBarry Smith 
2931baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W));
294c2fc9fa9SBarry Smith     */
295c2fc9fa9SBarry Smith     /* Compute a (scaled) negative update in the line search routine:
296c2fc9fa9SBarry Smith          Y <- X - lambda*Y
297c2fc9fa9SBarry Smith        and evaluate G = function(Y) (depends on the line search).
298c2fc9fa9SBarry Smith     */
2999566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, snes->vec_sol_update));
3009371c9d4SSatish Balay     ynorm = 1;
3019371c9d4SSatish Balay     gnorm = vi->phinorm;
3029566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y));
3039566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
3049566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm));
3059566063dSJacob 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));
306c2fc9fa9SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
307c2fc9fa9SBarry Smith     if (snes->domainerror) {
308c2fc9fa9SBarry Smith       snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
30922c6f798SBarry Smith       sdm->ops->computefunction = vi->computeuserfunction;
3103ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
311c2fc9fa9SBarry Smith     }
312422a814eSBarry Smith     if (lssucceed) {
313c2fc9fa9SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
314c2fc9fa9SBarry Smith         PetscBool ismin;
315c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3169566063dSJacob Faibussowitsch         PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin));
317c2fc9fa9SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
318c2fc9fa9SBarry Smith         break;
319c2fc9fa9SBarry Smith       }
320c2fc9fa9SBarry Smith     }
321c2fc9fa9SBarry Smith     /* Update function and solution vectors */
322c2fc9fa9SBarry Smith     vi->phinorm = gnorm;
323c2fc9fa9SBarry Smith     vi->merit   = 0.5 * vi->phinorm * vi->phinorm;
324c2fc9fa9SBarry Smith     /* Monitor convergence */
3259566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
326c2fc9fa9SBarry Smith     snes->iter  = i + 1;
327c2fc9fa9SBarry Smith     snes->norm  = vi->phinorm;
328c1e67a49SFande Kong     snes->xnorm = xnorm;
329c1e67a49SFande Kong     snes->ynorm = ynorm;
3309566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3319566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
332c2fc9fa9SBarry Smith     /* Test for convergence, xnorm = || X || */
3339566063dSJacob Faibussowitsch     if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
3342d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, vi->phinorm));
3352d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
336c2fc9fa9SBarry Smith     if (snes->reason) break;
337c2fc9fa9SBarry Smith   }
33822c6f798SBarry Smith   sdm->ops->computefunction = vi->computeuserfunction;
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
340c2fc9fa9SBarry Smith }
341c2fc9fa9SBarry Smith 
342d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes)
343d71ae5a4SJacob Faibussowitsch {
344f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
345c2fc9fa9SBarry Smith 
346c2fc9fa9SBarry Smith   PetscFunctionBegin;
3479566063dSJacob Faibussowitsch   PetscCall(SNESSetUp_VI(snes));
3489566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi));
3499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->phi));
3509566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->Da));
3519566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->Db));
3529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->z));
3539566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->t));
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355c2fc9fa9SBarry Smith }
356f6dfbefdSBarry Smith 
357d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes)
358d71ae5a4SJacob Faibussowitsch {
359f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
360c2fc9fa9SBarry Smith 
361c2fc9fa9SBarry Smith   PetscFunctionBegin;
3629566063dSJacob Faibussowitsch   PetscCall(SNESReset_VI(snes));
3639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->dpsi));
3649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->phi));
3659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->Da));
3669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->Db));
3679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->z));
3689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->t));
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
370c2fc9fa9SBarry Smith }
371c2fc9fa9SBarry Smith 
372d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject)
373d71ae5a4SJacob Faibussowitsch {
374c2fc9fa9SBarry Smith   PetscFunctionBegin;
375dbbe0bcdSBarry Smith   PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject));
376d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options");
377d0609cedSBarry Smith   PetscOptionsHeadEnd();
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
379c2fc9fa9SBarry Smith }
380c2fc9fa9SBarry Smith 
381c2fc9fa9SBarry Smith /*MC
382f450aa47SBarry Smith       SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method
383c2fc9fa9SBarry Smith 
384f6dfbefdSBarry Smith    Options Database Keys:
385d5f1b7e6SEd Bueler +   -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method
38661589011SJed Brown -   -snes_vi_monitor - prints the number of active constraints at each iteration.
38761589011SJed Brown 
388c2fc9fa9SBarry Smith    Level: beginner
389c2fc9fa9SBarry Smith 
390b80f3ac1SShri Abhyankar    References:
391606c0280SSatish Balay +  * -  T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth
392b80f3ac1SShri Abhyankar      algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001).
393606c0280SSatish Balay -  * -  T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale
394d5f1b7e6SEd Bueler      Applications, Optimization Methods and Software, 21 (2006).
395b80f3ac1SShri Abhyankar 
396f6dfbefdSBarry Smith    Notes:
397f6dfbefdSBarry Smith    This family of algorithm is much like an interior point method.
398c2fc9fa9SBarry Smith 
399f6dfbefdSBarry Smith    The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems
400f6dfbefdSBarry Smith 
401f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
402c2fc9fa9SBarry Smith M*/
403d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes)
404d71ae5a4SJacob Faibussowitsch {
405f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi;
406d8d34be6SBarry Smith   SNESLineSearch     linesearch;
407c2fc9fa9SBarry Smith 
408c2fc9fa9SBarry Smith   PetscFunctionBegin;
409f450aa47SBarry Smith   snes->ops->reset          = SNESReset_VINEWTONSSLS;
410f450aa47SBarry Smith   snes->ops->setup          = SNESSetUp_VINEWTONSSLS;
411f450aa47SBarry Smith   snes->ops->solve          = SNESSolve_VINEWTONSSLS;
412c2fc9fa9SBarry Smith   snes->ops->destroy        = SNESDestroy_VI;
413f450aa47SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS;
4140298fd71SBarry Smith   snes->ops->view           = NULL;
415c2fc9fa9SBarry Smith 
416c2fc9fa9SBarry Smith   snes->usesksp = PETSC_TRUE;
417efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
418c2fc9fa9SBarry Smith 
4199566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
420ec786807SJed Brown   if (!((PetscObject)linesearch)->type_name) {
4219566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
4229566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0));
423ec786807SJed Brown   }
424d8d34be6SBarry Smith 
4254fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
4264fc747eaSLawrence Mitchell 
4274dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&vi));
428c2fc9fa9SBarry Smith   snes->data = (void *)vi;
429c2fc9fa9SBarry Smith 
4309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI));
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI));
4323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
433c2fc9fa9SBarry Smith }
434