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