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