xref: /petsc/src/snes/impls/vi/ss/viss.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c2fc9fa9SBarry Smith 
2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/
3c2fc9fa9SBarry Smith 
4c2fc9fa9SBarry Smith /*
5c2fc9fa9SBarry Smith   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
6c2fc9fa9SBarry Smith 
7c2fc9fa9SBarry Smith   Input Parameter:
8c2fc9fa9SBarry Smith . phi - the semismooth function
9c2fc9fa9SBarry Smith 
10c2fc9fa9SBarry Smith   Output Parameter:
11c2fc9fa9SBarry Smith . merit - the merit function
12c2fc9fa9SBarry Smith . phinorm - ||phi||
13c2fc9fa9SBarry Smith 
14c2fc9fa9SBarry Smith   Notes:
15c2fc9fa9SBarry Smith   The merit function for the mixed complementarity problem is defined as
16c2fc9fa9SBarry Smith      merit = 0.5*phi^T*phi
17c2fc9fa9SBarry Smith */
18*9371c9d4SSatish Balay static 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 
23c2fc9fa9SBarry Smith   *merit = 0.5 * (*phinorm) * (*phinorm);
24c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
25c2fc9fa9SBarry Smith }
26c2fc9fa9SBarry Smith 
27*9371c9d4SSatish Balay static inline PetscScalar Phi(PetscScalar a, PetscScalar b) {
28c2fc9fa9SBarry Smith   return a + b - PetscSqrtScalar(a * a + b * b);
29c2fc9fa9SBarry Smith }
30c2fc9fa9SBarry Smith 
31*9371c9d4SSatish Balay static inline PetscScalar DPhi(PetscScalar a, PetscScalar b) {
32c2fc9fa9SBarry Smith   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b);
33c2fc9fa9SBarry Smith   else return .5;
34c2fc9fa9SBarry Smith }
35c2fc9fa9SBarry Smith 
36c2fc9fa9SBarry Smith /*
37c2fc9fa9SBarry Smith    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
38c2fc9fa9SBarry Smith 
39c2fc9fa9SBarry Smith    Input Parameters:
40c2fc9fa9SBarry Smith .  snes - the SNES context
41d5f1b7e6SEd Bueler .  X - current iterate
42c2fc9fa9SBarry Smith .  functx - user defined function context
43c2fc9fa9SBarry Smith 
44c2fc9fa9SBarry Smith    Output Parameters:
45c2fc9fa9SBarry Smith .  phi - Semismooth function
46c2fc9fa9SBarry Smith 
47c2fc9fa9SBarry Smith */
48*9371c9d4SSatish Balay static PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx) {
49f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
50c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu, F = snes->vec_func;
515edff71fSBarry Smith   PetscScalar       *phi_arr, *f_arr, *l, *u;
525edff71fSBarry Smith   const PetscScalar *x_arr;
53c2fc9fa9SBarry Smith   PetscInt           i, nlocal;
54c2fc9fa9SBarry Smith 
55c2fc9fa9SBarry Smith   PetscFunctionBegin;
569566063dSJacob Faibussowitsch   PetscCall((*vi->computeuserfunction)(snes, X, F, functx));
579566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x_arr));
599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f_arr));
609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xl, &l));
619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xu, &u));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(phi, &phi_arr));
63c2fc9fa9SBarry Smith 
64c2fc9fa9SBarry Smith   for (i = 0; i < nlocal; i++) {
65e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
66c2fc9fa9SBarry Smith       phi_arr[i] = f_arr[i];
67e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
68c2fc9fa9SBarry Smith       phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]);
69e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
70c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]);
71c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) {
72c2fc9fa9SBarry Smith       phi_arr[i] = l[i] - x_arr[i];
73c2fc9fa9SBarry Smith     } else { /* both bounds on variable */
74c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i]));
75c2fc9fa9SBarry Smith     }
76c2fc9fa9SBarry Smith   }
77c2fc9fa9SBarry Smith 
789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x_arr));
799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f_arr));
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xl, &l));
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xu, &u));
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(phi, &phi_arr));
83c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
84c2fc9fa9SBarry Smith }
85c2fc9fa9SBarry Smith 
86c2fc9fa9SBarry Smith /*
87c2fc9fa9SBarry Smith    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
88c2fc9fa9SBarry Smith                                           the semismooth jacobian.
89c2fc9fa9SBarry Smith */
90*9371c9d4SSatish Balay PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db) {
91c2fc9fa9SBarry Smith   PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2;
92c2fc9fa9SBarry Smith   PetscInt     i, nlocal;
93c2fc9fa9SBarry Smith 
94c2fc9fa9SBarry Smith   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
979566063dSJacob Faibussowitsch   PetscCall(VecGetArray(snes->xl, &l));
989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(snes->xu, &u));
999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Da, &da));
1009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Db, &db));
1019566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
102c2fc9fa9SBarry Smith 
103c2fc9fa9SBarry Smith   for (i = 0; i < nlocal; i++) {
104e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
105c2fc9fa9SBarry Smith       da[i] = 0;
106c2fc9fa9SBarry Smith       db[i] = 1;
107e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
108c2fc9fa9SBarry Smith       da[i] = DPhi(u[i] - x[i], -f[i]);
109c2fc9fa9SBarry Smith       db[i] = DPhi(-f[i], u[i] - x[i]);
110e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
111c2fc9fa9SBarry Smith       da[i] = DPhi(x[i] - l[i], f[i]);
112c2fc9fa9SBarry Smith       db[i] = DPhi(f[i], x[i] - l[i]);
113c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) { /* fixed variable */
114c2fc9fa9SBarry Smith       da[i] = 1;
115c2fc9fa9SBarry Smith       db[i] = 0;
116c2fc9fa9SBarry Smith     } else { /* upper and lower bounds on variable */
117c2fc9fa9SBarry Smith       da1   = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
118c2fc9fa9SBarry Smith       db1   = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]);
119c2fc9fa9SBarry Smith       da2   = DPhi(u[i] - x[i], -f[i]);
120c2fc9fa9SBarry Smith       db2   = DPhi(-f[i], u[i] - x[i]);
121c2fc9fa9SBarry Smith       da[i] = da1 + db1 * da2;
122c2fc9fa9SBarry Smith       db[i] = db1 * db2;
123c2fc9fa9SBarry Smith     }
124c2fc9fa9SBarry Smith   }
125c2fc9fa9SBarry Smith 
1269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(snes->xl, &l));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(snes->xu, &u));
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Da, &da));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Db, &db));
132c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
133c2fc9fa9SBarry Smith }
134c2fc9fa9SBarry Smith 
135c2fc9fa9SBarry Smith /*
136c2fc9fa9SBarry 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.
137c2fc9fa9SBarry Smith 
138c2fc9fa9SBarry Smith    Input Parameters:
139c2fc9fa9SBarry Smith .  Da       - Diagonal shift vector for the semismooth jacobian.
140c2fc9fa9SBarry Smith .  Db       - Row scaling vector for the semismooth jacobian.
141c2fc9fa9SBarry Smith 
142c2fc9fa9SBarry Smith    Output Parameters:
143c2fc9fa9SBarry Smith .  jac      - semismooth jacobian
144c2fc9fa9SBarry Smith .  jac_pre  - optional preconditioning matrix
145c2fc9fa9SBarry Smith 
146c2fc9fa9SBarry Smith    Notes:
147c2fc9fa9SBarry Smith    The semismooth jacobian matrix is given by
148c2fc9fa9SBarry Smith    jac = Da + Db*jacfun
149c2fc9fa9SBarry Smith    where Db is the row scaling matrix stored as a vector,
150c2fc9fa9SBarry Smith          Da is the diagonal perturbation matrix stored as a vector
151c2fc9fa9SBarry Smith    and   jacfun is the jacobian of the original nonlinear function.
152c2fc9fa9SBarry Smith */
153*9371c9d4SSatish Balay PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db) {
154c2fc9fa9SBarry Smith   /* Do row scaling  and add diagonal perturbation */
155362febeeSStefano Zampini   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(jac, Db, NULL));
1579566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES));
158c2fc9fa9SBarry Smith   if (jac != jac_pre) { /* If jac and jac_pre are different */
1599566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(jac_pre, Db, NULL));
1609566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES));
161c2fc9fa9SBarry Smith   }
162c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
163c2fc9fa9SBarry Smith }
164c2fc9fa9SBarry Smith 
165c2fc9fa9SBarry Smith /*
166c2fc9fa9SBarry Smith    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
167c2fc9fa9SBarry Smith 
168c2fc9fa9SBarry Smith    Input Parameters:
169c2fc9fa9SBarry Smith    phi - semismooth function.
170c2fc9fa9SBarry Smith    H   - semismooth jacobian
171c2fc9fa9SBarry Smith 
172c2fc9fa9SBarry Smith    Output Parameters:
173c2fc9fa9SBarry Smith    dpsi - merit function gradient
174c2fc9fa9SBarry Smith 
175c2fc9fa9SBarry Smith    Notes:
176c2fc9fa9SBarry Smith   The merit function gradient is computed as follows
177c2fc9fa9SBarry Smith         dpsi = H^T*phi
178c2fc9fa9SBarry Smith */
179*9371c9d4SSatish Balay PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) {
180c2fc9fa9SBarry Smith   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(H, phi, dpsi));
182c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
183c2fc9fa9SBarry Smith }
184c2fc9fa9SBarry Smith 
185c2fc9fa9SBarry Smith /*
186f450aa47SBarry Smith    SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton
187c2fc9fa9SBarry Smith    method using a line search.
188c2fc9fa9SBarry Smith 
189c2fc9fa9SBarry Smith    Input Parameters:
190c2fc9fa9SBarry Smith .  snes - the SNES context
191c2fc9fa9SBarry Smith 
192c2fc9fa9SBarry Smith    Application Interface Routine: SNESSolve()
193c2fc9fa9SBarry Smith 
194c2fc9fa9SBarry Smith    Notes:
195c2fc9fa9SBarry Smith    This implements essentially a semismooth Newton method with a
196d5f1b7e6SEd Bueler    line search. The default line search does not do any line search
197d5f1b7e6SEd Bueler    but rather takes a full Newton step.
198016d19f9SBarry Smith 
19904d7464bSBarry 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
200016d19f9SBarry Smith 
201c2fc9fa9SBarry Smith */
202*9371c9d4SSatish Balay PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) {
203f450aa47SBarry Smith   SNES_VINEWTONSSLS   *vi = (SNES_VINEWTONSSLS *)snes->data;
204c2fc9fa9SBarry Smith   PetscInt             maxits, i, lits;
205422a814eSBarry Smith   SNESLineSearchReason lssucceed;
206c2fc9fa9SBarry Smith   PetscReal            gnorm, xnorm = 0, ynorm;
2079bd66eb0SPeter Brune   Vec                  Y, X, F;
208c2fc9fa9SBarry Smith   KSPConvergedReason   kspreason;
2096cab3a1bSJed Brown   DM                   dm;
210942e3340SBarry Smith   DMSNES               sdm;
211c2fc9fa9SBarry Smith 
212c2fc9fa9SBarry Smith   PetscFunctionBegin;
2139566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
2149566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
2151aa26658SKarl Rupp 
21622c6f798SBarry Smith   vi->computeuserfunction   = sdm->ops->computefunction;
21722c6f798SBarry Smith   sdm->ops->computefunction = SNESVIComputeFunction;
218c2fc9fa9SBarry Smith 
219c2fc9fa9SBarry Smith   snes->numFailures            = 0;
220c2fc9fa9SBarry Smith   snes->numLinearSolveFailures = 0;
221c2fc9fa9SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
222c2fc9fa9SBarry Smith 
223c2fc9fa9SBarry Smith   maxits = snes->max_its;  /* maximum number of iterations */
224c2fc9fa9SBarry Smith   X      = snes->vec_sol;  /* solution vector */
225c2fc9fa9SBarry Smith   F      = snes->vec_func; /* residual vector */
226c2fc9fa9SBarry Smith   Y      = snes->work[0];  /* work vectors */
227c2fc9fa9SBarry Smith 
2289566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
229c2fc9fa9SBarry Smith   snes->iter = 0;
230c2fc9fa9SBarry Smith   snes->norm = 0.0;
2319566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
232c2fc9fa9SBarry Smith 
2339566063dSJacob Faibussowitsch   PetscCall(SNESVIProjectOntoBounds(snes, X));
2349566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, vi->phi));
235c2fc9fa9SBarry Smith   if (snes->domainerror) {
236c2fc9fa9SBarry Smith     snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
23722c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
238c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
239c2fc9fa9SBarry Smith   }
240c2fc9fa9SBarry Smith   /* Compute Merit function */
2419566063dSJacob Faibussowitsch   PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm));
242c2fc9fa9SBarry Smith 
2439566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x||  */
2449566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
245422a814eSBarry Smith   SNESCheckFunctionNorm(snes, vi->merit);
246c2fc9fa9SBarry Smith 
2479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
248c2fc9fa9SBarry Smith   snes->norm = vi->phinorm;
2499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2509566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0));
2519566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, vi->phinorm));
252c2fc9fa9SBarry Smith 
253c2fc9fa9SBarry Smith   /* test convergence */
254dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, vi->phinorm, &snes->reason, snes->cnvP);
255c2fc9fa9SBarry Smith   if (snes->reason) {
25622c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
257c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
258c2fc9fa9SBarry Smith   }
259c2fc9fa9SBarry Smith 
260c2fc9fa9SBarry Smith   for (i = 0; i < maxits; i++) {
261c2fc9fa9SBarry Smith     /* Call general purpose update function */
262dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
263c2fc9fa9SBarry Smith 
264c2fc9fa9SBarry Smith     /* Solve J Y = Phi, where J is the semismooth jacobian */
265b9600128SPeter Brune 
266b9600128SPeter Brune     /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/
26722c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
2689566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
26907b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
27022c6f798SBarry Smith     sdm->ops->computefunction = SNESVIComputeFunction;
271b9600128SPeter Brune 
272c2fc9fa9SBarry Smith     /* Get the diagonal shift and row scaling vectors */
2739566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db));
274c2fc9fa9SBarry Smith     /* Compute the semismooth jacobian */
2759566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db));
276c2fc9fa9SBarry Smith     /* Compute the merit function gradient */
2779566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi));
2789566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2799566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, vi->phi, Y));
2809566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason));
281c2fc9fa9SBarry Smith 
282c2fc9fa9SBarry Smith     if (kspreason < 0) {
283c2fc9fa9SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
28463a3b9bcSJacob 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));
285c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
286c2fc9fa9SBarry Smith         break;
287c2fc9fa9SBarry Smith       }
288c2fc9fa9SBarry Smith     }
2899566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
290c2fc9fa9SBarry Smith     snes->linear_its += lits;
29163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
292c2fc9fa9SBarry Smith     /*
2936b2b7091SBarry Smith     if (snes->ops->precheck) {
294c2fc9fa9SBarry Smith       PetscBool changed_y = PETSC_FALSE;
295dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y);
296c2fc9fa9SBarry Smith     }
297c2fc9fa9SBarry Smith 
2981baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W));
299c2fc9fa9SBarry Smith     */
300c2fc9fa9SBarry Smith     /* Compute a (scaled) negative update in the line search routine:
301c2fc9fa9SBarry Smith          Y <- X - lambda*Y
302c2fc9fa9SBarry Smith        and evaluate G = function(Y) (depends on the line search).
303c2fc9fa9SBarry Smith     */
3049566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, snes->vec_sol_update));
305*9371c9d4SSatish Balay     ynorm = 1;
306*9371c9d4SSatish Balay     gnorm = vi->phinorm;
3079566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y));
3089566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
3099566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm));
3109566063dSJacob 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));
311c2fc9fa9SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
312c2fc9fa9SBarry Smith     if (snes->domainerror) {
313c2fc9fa9SBarry Smith       snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
31422c6f798SBarry Smith       sdm->ops->computefunction = vi->computeuserfunction;
315c2fc9fa9SBarry Smith       PetscFunctionReturn(0);
316c2fc9fa9SBarry Smith     }
317422a814eSBarry Smith     if (lssucceed) {
318c2fc9fa9SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
319c2fc9fa9SBarry Smith         PetscBool ismin;
320c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3219566063dSJacob Faibussowitsch         PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin));
322c2fc9fa9SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
323c2fc9fa9SBarry Smith         break;
324c2fc9fa9SBarry Smith       }
325c2fc9fa9SBarry Smith     }
326c2fc9fa9SBarry Smith     /* Update function and solution vectors */
327c2fc9fa9SBarry Smith     vi->phinorm = gnorm;
328c2fc9fa9SBarry Smith     vi->merit   = 0.5 * vi->phinorm * vi->phinorm;
329c2fc9fa9SBarry Smith     /* Monitor convergence */
3309566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
331c2fc9fa9SBarry Smith     snes->iter  = i + 1;
332c2fc9fa9SBarry Smith     snes->norm  = vi->phinorm;
333c1e67a49SFande Kong     snes->xnorm = xnorm;
334c1e67a49SFande Kong     snes->ynorm = ynorm;
3359566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3369566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
3379566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
338c2fc9fa9SBarry Smith     /* Test for convergence, xnorm = || X || */
3399566063dSJacob Faibussowitsch     if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
340dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, vi->phinorm, &snes->reason, snes->cnvP);
341c2fc9fa9SBarry Smith     if (snes->reason) break;
342c2fc9fa9SBarry Smith   }
343c2fc9fa9SBarry Smith   if (i == maxits) {
34463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
345c2fc9fa9SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
346c2fc9fa9SBarry Smith   }
34722c6f798SBarry Smith   sdm->ops->computefunction = vi->computeuserfunction;
348c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
349c2fc9fa9SBarry Smith }
350c2fc9fa9SBarry Smith 
351c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
352c2fc9fa9SBarry Smith /*
353f450aa47SBarry Smith    SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use
354c2fc9fa9SBarry Smith    of the SNES nonlinear solver.
355c2fc9fa9SBarry Smith 
356c2fc9fa9SBarry Smith    Input Parameter:
357c2fc9fa9SBarry Smith .  snes - the SNES context
358c2fc9fa9SBarry Smith 
359c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetUp()
360c2fc9fa9SBarry Smith 
361c2fc9fa9SBarry Smith    Notes:
362c2fc9fa9SBarry Smith    For basic use of the SNES solvers, the user need not explicitly call
363c2fc9fa9SBarry Smith    SNESSetUp(), since these actions will automatically occur during
364c2fc9fa9SBarry Smith    the call to SNESSolve().
365c2fc9fa9SBarry Smith  */
366*9371c9d4SSatish Balay PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) {
367f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
368c2fc9fa9SBarry Smith 
369c2fc9fa9SBarry Smith   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   PetscCall(SNESSetUp_VI(snes));
3719566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi));
3729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->phi));
3739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->Da));
3749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->Db));
3759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->z));
3769566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(snes->vec_sol, &vi->t));
377c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
378c2fc9fa9SBarry Smith }
379c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
380*9371c9d4SSatish Balay PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) {
381f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
382c2fc9fa9SBarry Smith 
383c2fc9fa9SBarry Smith   PetscFunctionBegin;
3849566063dSJacob Faibussowitsch   PetscCall(SNESReset_VI(snes));
3859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->dpsi));
3869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->phi));
3879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->Da));
3889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->Db));
3899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->z));
3909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vi->t));
391c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
392c2fc9fa9SBarry Smith }
393c2fc9fa9SBarry Smith 
394c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
395c2fc9fa9SBarry Smith /*
396f450aa47SBarry Smith    SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method.
397c2fc9fa9SBarry Smith 
398c2fc9fa9SBarry Smith    Input Parameter:
399c2fc9fa9SBarry Smith .  snes - the SNES context
400c2fc9fa9SBarry Smith 
401c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetFromOptions()
402c2fc9fa9SBarry Smith */
403*9371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject) {
404c2fc9fa9SBarry Smith   PetscFunctionBegin;
405dbbe0bcdSBarry Smith   PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject));
406d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options");
407d0609cedSBarry Smith   PetscOptionsHeadEnd();
408c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
409c2fc9fa9SBarry Smith }
410c2fc9fa9SBarry Smith 
411c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
412c2fc9fa9SBarry Smith /*MC
413f450aa47SBarry Smith       SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method
414c2fc9fa9SBarry Smith 
41561589011SJed Brown    Options Database:
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 
427c2e3fba1SPatrick Sanan .seealso: `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
428c2fc9fa9SBarry Smith 
429c2fc9fa9SBarry Smith M*/
430*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) {
431f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi;
432d8d34be6SBarry Smith   SNESLineSearch     linesearch;
433c2fc9fa9SBarry Smith 
434c2fc9fa9SBarry Smith   PetscFunctionBegin;
435f450aa47SBarry Smith   snes->ops->reset          = SNESReset_VINEWTONSSLS;
436f450aa47SBarry Smith   snes->ops->setup          = SNESSetUp_VINEWTONSSLS;
437f450aa47SBarry Smith   snes->ops->solve          = SNESSolve_VINEWTONSSLS;
438c2fc9fa9SBarry Smith   snes->ops->destroy        = SNESDestroy_VI;
439f450aa47SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS;
4400298fd71SBarry Smith   snes->ops->view           = NULL;
441c2fc9fa9SBarry Smith 
442c2fc9fa9SBarry Smith   snes->usesksp = PETSC_TRUE;
443efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
444c2fc9fa9SBarry Smith 
4459566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
446ec786807SJed Brown   if (!((PetscObject)linesearch)->type_name) {
4479566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
4489566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0));
449ec786807SJed Brown   }
450d8d34be6SBarry Smith 
4514fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
4524fc747eaSLawrence Mitchell 
4539566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &vi));
454c2fc9fa9SBarry Smith   snes->data = (void *)vi;
455c2fc9fa9SBarry Smith 
4569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI));
4579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI));
458c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
459c2fc9fa9SBarry Smith }
460