xref: /petsc/src/snes/impls/vi/ss/viss.c (revision bcee047adeeb73090d7e36cc71e39fc287cdbb97)
1 
2 #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/
3 
4 /*@
5   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
6 
7    Input Parameter:
8 .  phi - the `Vec` holding the evaluation of the semismooth function
9 
10    Output Parameters:
11 +  merit - the merit function 1/2 ||phi||^2
12 -  phinorm -  the two-norm of the vector, ||phi||
13 
14    Level: developer
15 
16 .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeFunction()`
17 @*/
18 PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm)
19 {
20   PetscFunctionBegin;
21   PetscCall(VecNormBegin(phi, NORM_2, phinorm));
22   PetscCall(VecNormEnd(phi, NORM_2, phinorm));
23   *merit = 0.5 * (*phinorm) * (*phinorm);
24   PetscFunctionReturn(PETSC_SUCCESS);
25 }
26 
27 static inline PetscScalar Phi(PetscScalar a, PetscScalar b)
28 {
29   return a + b - PetscSqrtScalar(a * a + b * b);
30 }
31 
32 static inline PetscScalar DPhi(PetscScalar a, PetscScalar b)
33 {
34   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b);
35   else return .5;
36 }
37 
38 /*@
39    SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear
40    equations in semismooth form.
41 
42    Input Parameters:
43 +  snes - the SNES context
44 .  X - current iterate
45 -  functx - user defined function context
46 
47    Output Parameter:
48 .  phi - the evaluation of Semismooth function at X
49 
50    Level: developer
51 
52 .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()`
53 @*/
54 PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx)
55 {
56   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
57   Vec                Xl = snes->xl, Xu = snes->xu, F = snes->vec_func;
58   PetscScalar       *phi_arr, *f_arr, *l, *u;
59   const PetscScalar *x_arr;
60   PetscInt           i, nlocal;
61 
62   PetscFunctionBegin;
63   PetscCall((*vi->computeuserfunction)(snes, X, F, functx));
64   PetscCall(VecGetLocalSize(X, &nlocal));
65   PetscCall(VecGetArrayRead(X, &x_arr));
66   PetscCall(VecGetArray(F, &f_arr));
67   PetscCall(VecGetArray(Xl, &l));
68   PetscCall(VecGetArray(Xu, &u));
69   PetscCall(VecGetArray(phi, &phi_arr));
70 
71   for (i = 0; i < nlocal; i++) {
72     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
73       phi_arr[i] = f_arr[i];
74     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
75       phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]);
76     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
77       phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]);
78     } else if (l[i] == u[i]) {
79       phi_arr[i] = l[i] - x_arr[i];
80     } else { /* both bounds on variable */
81       phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i]));
82     }
83   }
84 
85   PetscCall(VecRestoreArrayRead(X, &x_arr));
86   PetscCall(VecRestoreArray(F, &f_arr));
87   PetscCall(VecRestoreArray(Xl, &l));
88   PetscCall(VecRestoreArray(Xu, &u));
89   PetscCall(VecRestoreArray(phi, &phi_arr));
90   PetscFunctionReturn(PETSC_SUCCESS);
91 }
92 
93 /*
94    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
95                                           the semismooth jacobian.
96 */
97 PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db)
98 {
99   PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2;
100   PetscInt     i, nlocal;
101 
102   PetscFunctionBegin;
103   PetscCall(VecGetArray(X, &x));
104   PetscCall(VecGetArray(F, &f));
105   PetscCall(VecGetArray(snes->xl, &l));
106   PetscCall(VecGetArray(snes->xu, &u));
107   PetscCall(VecGetArray(Da, &da));
108   PetscCall(VecGetArray(Db, &db));
109   PetscCall(VecGetLocalSize(X, &nlocal));
110 
111   for (i = 0; i < nlocal; i++) {
112     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
113       da[i] = 0;
114       db[i] = 1;
115     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */
116       da[i] = DPhi(u[i] - x[i], -f[i]);
117       db[i] = DPhi(-f[i], u[i] - x[i]);
118     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */
119       da[i] = DPhi(x[i] - l[i], f[i]);
120       db[i] = DPhi(f[i], x[i] - l[i]);
121     } else if (l[i] == u[i]) { /* fixed variable */
122       da[i] = 1;
123       db[i] = 0;
124     } else { /* upper and lower bounds on variable */
125       da1   = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
126       db1   = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]);
127       da2   = DPhi(u[i] - x[i], -f[i]);
128       db2   = DPhi(-f[i], u[i] - x[i]);
129       da[i] = da1 + db1 * da2;
130       db[i] = db1 * db2;
131     }
132   }
133 
134   PetscCall(VecRestoreArray(X, &x));
135   PetscCall(VecRestoreArray(F, &f));
136   PetscCall(VecRestoreArray(snes->xl, &l));
137   PetscCall(VecRestoreArray(snes->xu, &u));
138   PetscCall(VecRestoreArray(Da, &da));
139   PetscCall(VecRestoreArray(Db, &db));
140   PetscFunctionReturn(PETSC_SUCCESS);
141 }
142 
143 /*
144    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.
145 
146    Input Parameters:
147 .  Da       - Diagonal shift vector for the semismooth jacobian.
148 .  Db       - Row scaling vector for the semismooth jacobian.
149 
150    Output Parameters:
151 .  jac      - semismooth jacobian
152 .  jac_pre  - optional preconditioning matrix
153 
154    Note:
155    The semismooth jacobian matrix is given by
156    jac = Da + Db*jacfun
157    where Db is the row scaling matrix stored as a vector,
158          Da is the diagonal perturbation matrix stored as a vector
159    and   jacfun is the jacobian of the original nonlinear function.
160 */
161 PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db)
162 {
163   /* Do row scaling  and add diagonal perturbation */
164   PetscFunctionBegin;
165   PetscCall(MatDiagonalScale(jac, Db, NULL));
166   PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES));
167   if (jac != jac_pre) { /* If jac and jac_pre are different */
168     PetscCall(MatDiagonalScale(jac_pre, Db, NULL));
169     PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES));
170   }
171   PetscFunctionReturn(PETSC_SUCCESS);
172 }
173 
174 /*
175    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
176 
177    Input Parameters:
178    phi - semismooth function.
179    H   - semismooth jacobian
180 
181    Output Parameter:
182    dpsi - merit function gradient
183 
184    Note:
185   The merit function gradient is computed as follows
186         dpsi = H^T*phi
187 */
188 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
189 {
190   PetscFunctionBegin;
191   PetscCall(MatMultTranspose(H, phi, dpsi));
192   PetscFunctionReturn(PETSC_SUCCESS);
193 }
194 
195 /*
196    SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton
197    method using a line search.
198 
199    Input Parameter:
200 .  snes - the SNES context
201 
202    Application Interface Routine: SNESSolve()
203 
204    Note:
205    This implements essentially a semismooth Newton method with a
206    line search. The default line search does not do any line search
207    but rather takes a full Newton step.
208 
209    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
210 
211 */
212 PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes)
213 {
214   SNES_VINEWTONSSLS   *vi = (SNES_VINEWTONSSLS *)snes->data;
215   PetscInt             maxits, i, lits;
216   SNESLineSearchReason lssucceed;
217   PetscReal            gnorm, xnorm = 0, ynorm;
218   Vec                  Y, X, F;
219   KSPConvergedReason   kspreason;
220   DM                   dm;
221   DMSNES               sdm;
222 
223   PetscFunctionBegin;
224   PetscCall(SNESGetDM(snes, &dm));
225   PetscCall(DMGetDMSNES(dm, &sdm));
226 
227   vi->computeuserfunction   = sdm->ops->computefunction;
228   sdm->ops->computefunction = SNESVIComputeFunction;
229 
230   snes->numFailures            = 0;
231   snes->numLinearSolveFailures = 0;
232   snes->reason                 = SNES_CONVERGED_ITERATING;
233 
234   maxits = snes->max_its;  /* maximum number of iterations */
235   X      = snes->vec_sol;  /* solution vector */
236   F      = snes->vec_func; /* residual vector */
237   Y      = snes->work[0];  /* work vectors */
238 
239   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
240   snes->iter = 0;
241   snes->norm = 0.0;
242   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
243 
244   PetscCall(SNESVIProjectOntoBounds(snes, X));
245   PetscCall(SNESComputeFunction(snes, X, vi->phi));
246   if (snes->domainerror) {
247     snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
248     sdm->ops->computefunction = vi->computeuserfunction;
249     PetscFunctionReturn(PETSC_SUCCESS);
250   }
251   /* Compute Merit function */
252   PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm));
253 
254   PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x||  */
255   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
256   SNESCheckFunctionNorm(snes, vi->merit);
257 
258   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
259   snes->norm = vi->phinorm;
260   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
261   PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0));
262 
263   /* test convergence */
264   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, vi->phinorm));
265   PetscCall(SNESMonitor(snes, 0, vi->phinorm));
266   if (snes->reason) {
267     sdm->ops->computefunction = vi->computeuserfunction;
268     PetscFunctionReturn(PETSC_SUCCESS);
269   }
270 
271   for (i = 0; i < maxits; i++) {
272     /* Call general purpose update function */
273     PetscTryTypeMethod(snes, update, snes->iter);
274 
275     /* Solve J Y = Phi, where J is the semismooth jacobian */
276 
277     /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/
278     sdm->ops->computefunction = vi->computeuserfunction;
279     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
280     SNESCheckJacobianDomainerror(snes);
281     sdm->ops->computefunction = SNESVIComputeFunction;
282 
283     /* Get the diagonal shift and row scaling vectors */
284     PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db));
285     /* Compute the semismooth jacobian */
286     PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db));
287     /* Compute the merit function gradient */
288     PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi));
289     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
290     PetscCall(KSPSolve(snes->ksp, vi->phi, Y));
291     PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason));
292 
293     if (kspreason < 0) {
294       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
295         PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures));
296         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
297         break;
298       }
299     }
300     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
301     snes->linear_its += lits;
302     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
303     /*
304     if (snes->ops->precheck) {
305       PetscBool changed_y = PETSC_FALSE;
306       PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y);
307     }
308 
309     if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W));
310     */
311     /* Compute a (scaled) negative update in the line search routine:
312          Y <- X - lambda*Y
313        and evaluate G = function(Y) (depends on the line search).
314     */
315     PetscCall(VecCopy(Y, snes->vec_sol_update));
316     ynorm = 1;
317     gnorm = vi->phinorm;
318     PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y));
319     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
320     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm));
321     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)vi->phinorm, (double)gnorm, (double)ynorm, (int)lssucceed));
322     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
323     if (snes->domainerror) {
324       snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
325       sdm->ops->computefunction = vi->computeuserfunction;
326       PetscFunctionReturn(PETSC_SUCCESS);
327     }
328     if (lssucceed) {
329       if (++snes->numFailures >= snes->maxFailures) {
330         PetscBool ismin;
331         snes->reason = SNES_DIVERGED_LINE_SEARCH;
332         PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin));
333         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
334         break;
335       }
336     }
337     /* Update function and solution vectors */
338     vi->phinorm = gnorm;
339     vi->merit   = 0.5 * vi->phinorm * vi->phinorm;
340     /* Monitor convergence */
341     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
342     snes->iter  = i + 1;
343     snes->norm  = vi->phinorm;
344     snes->xnorm = xnorm;
345     snes->ynorm = ynorm;
346     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
347     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
348     /* Test for convergence, xnorm = || X || */
349     if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
350     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, vi->phinorm));
351     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
352     if (snes->reason) break;
353   }
354   sdm->ops->computefunction = vi->computeuserfunction;
355   PetscFunctionReturn(PETSC_SUCCESS);
356 }
357 
358 /*
359    SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use
360    of the SNES nonlinear solver.
361 
362    Input Parameter:
363 .  snes - the SNES context
364 
365    Application Interface Routine: SNESSetUp()
366 
367    Note:
368    For basic use of the SNES solvers, the user need not explicitly call
369    SNESSetUp(), since these actions will automatically occur during
370    the call to SNESSolve().
371  */
372 PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes)
373 {
374   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
375 
376   PetscFunctionBegin;
377   PetscCall(SNESSetUp_VI(snes));
378   PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi));
379   PetscCall(VecDuplicate(snes->vec_sol, &vi->phi));
380   PetscCall(VecDuplicate(snes->vec_sol, &vi->Da));
381   PetscCall(VecDuplicate(snes->vec_sol, &vi->Db));
382   PetscCall(VecDuplicate(snes->vec_sol, &vi->z));
383   PetscCall(VecDuplicate(snes->vec_sol, &vi->t));
384   PetscFunctionReturn(PETSC_SUCCESS);
385 }
386 
387 PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes)
388 {
389   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data;
390 
391   PetscFunctionBegin;
392   PetscCall(SNESReset_VI(snes));
393   PetscCall(VecDestroy(&vi->dpsi));
394   PetscCall(VecDestroy(&vi->phi));
395   PetscCall(VecDestroy(&vi->Da));
396   PetscCall(VecDestroy(&vi->Db));
397   PetscCall(VecDestroy(&vi->z));
398   PetscCall(VecDestroy(&vi->t));
399   PetscFunctionReturn(PETSC_SUCCESS);
400 }
401 
402 /*
403    SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method.
404 
405    Input Parameter:
406 .  snes - the SNES context
407 
408    Application Interface Routine: SNESSetFromOptions()
409 */
410 static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject)
411 {
412   PetscFunctionBegin;
413   PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject));
414   PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options");
415   PetscOptionsHeadEnd();
416   PetscFunctionReturn(PETSC_SUCCESS);
417 }
418 
419 /*MC
420       SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method
421 
422    Options Database Keys:
423 +   -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method
424 -   -snes_vi_monitor - prints the number of active constraints at each iteration.
425 
426    Level: beginner
427 
428    References:
429 +  * -  T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth
430      algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001).
431 -  * -  T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale
432      Applications, Optimization Methods and Software, 21 (2006).
433 
434    Notes:
435    This family of algorithm is much like an interior point method.
436 
437    The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems
438 
439 .seealso: `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
440 M*/
441 PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes)
442 {
443   SNES_VINEWTONSSLS *vi;
444   SNESLineSearch     linesearch;
445 
446   PetscFunctionBegin;
447   snes->ops->reset          = SNESReset_VINEWTONSSLS;
448   snes->ops->setup          = SNESSetUp_VINEWTONSSLS;
449   snes->ops->solve          = SNESSolve_VINEWTONSSLS;
450   snes->ops->destroy        = SNESDestroy_VI;
451   snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS;
452   snes->ops->view           = NULL;
453 
454   snes->usesksp = PETSC_TRUE;
455   snes->usesnpc = PETSC_FALSE;
456 
457   PetscCall(SNESGetLineSearch(snes, &linesearch));
458   if (!((PetscObject)linesearch)->type_name) {
459     PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
460     PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0));
461   }
462 
463   snes->alwayscomputesfinalresidual = PETSC_FALSE;
464 
465   PetscCall(PetscNew(&vi));
466   snes->data = (void *)vi;
467 
468   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI));
469   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI));
470   PetscFunctionReturn(PETSC_SUCCESS);
471 }
472