xref: /petsc/src/tao/util/tao_util.c (revision 174dc0c8cee294b82b85e4dd3b331b29396264fc)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/
38370d7cdSHansol Suh #include <petscsys.h>
4a7e14dcfSSatish Balay 
5d71ae5a4SJacob Faibussowitsch static inline PetscReal Fischer(PetscReal a, PetscReal b)
6d71ae5a4SJacob Faibussowitsch {
7a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
8ad540459SPierre Jolivet   if (a + b <= 0) return PetscSqrtReal(a * a + b * b) - (a + b);
946bdf8c8SLisandro Dalcin   return -2.0 * a * b / (PetscSqrtReal(a * a + b * b) + (a + b));
10a7e14dcfSSatish Balay }
11a7e14dcfSSatish Balay 
12a7e14dcfSSatish Balay /*@
13a7e14dcfSSatish Balay   VecFischer - Evaluates the Fischer-Burmeister function for complementarity
14a7e14dcfSSatish Balay   problems.
15a7e14dcfSSatish Balay 
16c3339decSBarry Smith   Logically Collective
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay   Input Parameters:
19a7e14dcfSSatish Balay + X - current point
20a7e14dcfSSatish Balay . F - function evaluated at x
21a7e14dcfSSatish Balay . L - lower bounds
22a7e14dcfSSatish Balay - U - upper bounds
23a7e14dcfSSatish Balay 
24f899ff85SJose E. Roman   Output Parameter:
25a7e14dcfSSatish Balay . FB - The Fischer-Burmeister function vector
26a7e14dcfSSatish Balay 
272fe279fdSBarry Smith   Level: developer
282fe279fdSBarry Smith 
29a7e14dcfSSatish Balay   Notes:
30a7e14dcfSSatish Balay   The Fischer-Burmeister function is defined as
31*b44f4de4SBarry Smith 
32*b44f4de4SBarry Smith   $$
33*b44f4de4SBarry Smith   \phi(a,b) := \sqrt{a^2 + b^2} - a - b
34*b44f4de4SBarry Smith   $$
35*b44f4de4SBarry Smith 
36a7e14dcfSSatish Balay   and is used reformulate a complementarity problem as a semismooth
37a7e14dcfSSatish Balay   system of equations.
38a7e14dcfSSatish Balay 
393b242c63SJacob Faibussowitsch   The result of this function is done by cases\:
40a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i]
41a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i])
42a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i])
43a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u]))
44a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
45a7e14dcfSSatish Balay 
46a1cb98faSBarry Smith .seealso: `Vec`, `VecSFischer()`, `MatDFischer()`, `MatDSFischer()`
47a7e14dcfSSatish Balay @*/
48d71ae5a4SJacob Faibussowitsch PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB)
49d71ae5a4SJacob Faibussowitsch {
5046bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
5146bdf8c8SLisandro Dalcin   PetscScalar       *fb;
52a7e14dcfSSatish Balay   PetscReal          xval, fval, lval, uval;
53a7e14dcfSSatish Balay   PetscInt           low[5], high[5], n, i;
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay   PetscFunctionBegin;
56a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID, 1);
57a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID, 2);
5876be6f4fSStefano Zampini   if (L) PetscValidHeaderSpecific(L, VEC_CLASSID, 3);
5976be6f4fSStefano Zampini   if (U) PetscValidHeaderSpecific(U, VEC_CLASSID, 4);
60064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(FB, VEC_CLASSID, 5);
61a7e14dcfSSatish Balay 
6276be6f4fSStefano Zampini   if (!L && !U) {
6376be6f4fSStefano Zampini     PetscCall(VecAXPBY(FB, -1.0, 0.0, F));
643ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6576be6f4fSStefano Zampini   }
6676be6f4fSStefano Zampini 
679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, low, high));
689566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(F, low + 1, high + 1));
699566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(L, low + 2, high + 2));
709566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(U, low + 3, high + 3));
719566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4));
72a7e14dcfSSatish Balay 
73ad540459SPierre Jolivet   for (i = 1; i < 4; ++i) PetscCheck(low[0] == low[i] && high[0] == high[i], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Vectors must be identically loaded over processors");
74a7e14dcfSSatish Balay 
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(L, &l));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(FB, &fb));
80a7e14dcfSSatish Balay 
819566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
82a7e14dcfSSatish Balay 
83a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
8476be6f4fSStefano Zampini     xval = PetscRealPart(x[i]);
8576be6f4fSStefano Zampini     fval = PetscRealPart(f[i]);
8676be6f4fSStefano Zampini     lval = PetscRealPart(l[i]);
8776be6f4fSStefano Zampini     uval = PetscRealPart(u[i]);
88a7e14dcfSSatish Balay 
8976be6f4fSStefano Zampini     if (lval <= -PETSC_INFINITY && uval >= PETSC_INFINITY) {
90a7e14dcfSSatish Balay       fb[i] = -fval;
91e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
92a7e14dcfSSatish Balay       fb[i] = -Fischer(uval - xval, -fval);
93e270355aSBarry Smith     } else if (uval >= PETSC_INFINITY) {
94a7e14dcfSSatish Balay       fb[i] = Fischer(xval - lval, fval);
952d0e5244SBarry Smith     } else if (lval == uval) {
96a7e14dcfSSatish Balay       fb[i] = lval - xval;
972d0e5244SBarry Smith     } else {
98a7e14dcfSSatish Balay       fval  = Fischer(uval - xval, -fval);
99a7e14dcfSSatish Balay       fb[i] = Fischer(xval - lval, fval);
100a7e14dcfSSatish Balay     }
101a7e14dcfSSatish Balay   }
102a7e14dcfSSatish Balay 
1039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
1059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(L, &l));
1069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(FB, &fb));
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109a7e14dcfSSatish Balay }
110a7e14dcfSSatish Balay 
111d71ae5a4SJacob Faibussowitsch static inline PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c)
112d71ae5a4SJacob Faibussowitsch {
113a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
114ad540459SPierre Jolivet   if (a + b <= 0) return PetscSqrtReal(a * a + b * b + 2.0 * c * c) - (a + b);
1153f6ba705SLisandro Dalcin   return 2.0 * (c * c - a * b) / (PetscSqrtReal(a * a + b * b + 2.0 * c * c) + (a + b));
116a7e14dcfSSatish Balay }
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay /*@
119a7e14dcfSSatish Balay   VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for
120a7e14dcfSSatish Balay   complementarity problems.
121a7e14dcfSSatish Balay 
122c3339decSBarry Smith   Logically Collective
123a7e14dcfSSatish Balay 
124a7e14dcfSSatish Balay   Input Parameters:
125a7e14dcfSSatish Balay + X  - current point
126a7e14dcfSSatish Balay . F  - function evaluated at x
127a7e14dcfSSatish Balay . L  - lower bounds
128a7e14dcfSSatish Balay . U  - upper bounds
129a7e14dcfSSatish Balay - mu - smoothing parameter
130a7e14dcfSSatish Balay 
131f899ff85SJose E. Roman   Output Parameter:
132a7e14dcfSSatish Balay . FB - The Smoothed Fischer-Burmeister function vector
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay   Notes:
135a7e14dcfSSatish Balay   The Smoothed Fischer-Burmeister function is defined as
136a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b
137a7e14dcfSSatish Balay   and is used reformulate a complementarity problem as a semismooth
138a7e14dcfSSatish Balay   system of equations.
139a7e14dcfSSatish Balay 
1403b242c63SJacob Faibussowitsch   The result of this function is done by cases\:
141a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i] - 2*mu*x[i]
142a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i], mu)
143a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i], mu)
144a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu)
145a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay   Level: developer
148a7e14dcfSSatish Balay 
149a1cb98faSBarry Smith .seealso: `Vec`, `VecFischer()`, `MatDFischer()`, `MatDSFischer()`
150a7e14dcfSSatish Balay @*/
151d71ae5a4SJacob Faibussowitsch PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB)
152d71ae5a4SJacob Faibussowitsch {
15346bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
15446bdf8c8SLisandro Dalcin   PetscScalar       *fb;
155a7e14dcfSSatish Balay   PetscReal          xval, fval, lval, uval;
156a7e14dcfSSatish Balay   PetscInt           low[5], high[5], n, i;
157a7e14dcfSSatish Balay 
158a7e14dcfSSatish Balay   PetscFunctionBegin;
159a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID, 1);
160a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID, 2);
161a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID, 3);
162a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID, 4);
163a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID, 6);
164a7e14dcfSSatish Balay 
1659566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, low, high));
1669566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(F, low + 1, high + 1));
1679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(L, low + 2, high + 2));
1689566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(U, low + 3, high + 3));
1699566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4));
170a7e14dcfSSatish Balay 
171ad540459SPierre Jolivet   for (i = 1; i < 4; ++i) PetscCheck(low[0] == low[i] && high[0] == high[i], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Vectors must be identically loaded over processors");
172a7e14dcfSSatish Balay 
1739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
1759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(L, &l));
1769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(FB, &fb));
178a7e14dcfSSatish Balay 
1799566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
1829371c9d4SSatish Balay     xval = PetscRealPart(*x++);
1839371c9d4SSatish Balay     fval = PetscRealPart(*f++);
1849371c9d4SSatish Balay     lval = PetscRealPart(*l++);
1859371c9d4SSatish Balay     uval = PetscRealPart(*u++);
186a7e14dcfSSatish Balay 
187e270355aSBarry Smith     if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) {
188a7e14dcfSSatish Balay       (*fb++) = -fval - mu * xval;
189e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
190a7e14dcfSSatish Balay       (*fb++) = -SFischer(uval - xval, -fval, mu);
191e270355aSBarry Smith     } else if (uval >= PETSC_INFINITY) {
192a7e14dcfSSatish Balay       (*fb++) = SFischer(xval - lval, fval, mu);
1932d0e5244SBarry Smith     } else if (lval == uval) {
194a7e14dcfSSatish Balay       (*fb++) = lval - xval;
1952d0e5244SBarry Smith     } else {
196a7e14dcfSSatish Balay       fval    = SFischer(uval - xval, -fval, mu);
197a7e14dcfSSatish Balay       (*fb++) = SFischer(xval - lval, fval, mu);
198a7e14dcfSSatish Balay     }
199a7e14dcfSSatish Balay   }
2009371c9d4SSatish Balay   x -= n;
2019371c9d4SSatish Balay   f -= n;
2029371c9d4SSatish Balay   l -= n;
2039371c9d4SSatish Balay   u -= n;
2049371c9d4SSatish Balay   fb -= n;
205a7e14dcfSSatish Balay 
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(L, &l));
2099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(FB, &fb));
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212a7e14dcfSSatish Balay }
213a7e14dcfSSatish Balay 
214d71ae5a4SJacob Faibussowitsch static inline PetscReal fischnorm(PetscReal a, PetscReal b)
215d71ae5a4SJacob Faibussowitsch {
216658c1fc4SLisandro Dalcin   return PetscSqrtReal(a * a + b * b);
217a7e14dcfSSatish Balay }
218a7e14dcfSSatish Balay 
219d71ae5a4SJacob Faibussowitsch static inline PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c)
220d71ae5a4SJacob Faibussowitsch {
221658c1fc4SLisandro Dalcin   return PetscSqrtReal(a * a + b * b + 2.0 * c * c);
222a7e14dcfSSatish Balay }
223a7e14dcfSSatish Balay 
224a7e14dcfSSatish Balay /*@
225235fd6e6SBarry Smith   MatDFischer - Calculates an element of the B-subdifferential of the
226a7e14dcfSSatish Balay   Fischer-Burmeister function for complementarity problems.
227a7e14dcfSSatish Balay 
228c3339decSBarry Smith   Collective
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay   Input Parameters:
2313b242c63SJacob Faibussowitsch + jac - the jacobian of `f` at `X`
232a7e14dcfSSatish Balay . X   - current point
2333b242c63SJacob Faibussowitsch . Con - constraints function evaluated at `X`
234a7e14dcfSSatish Balay . XL  - lower bounds
235a7e14dcfSSatish Balay . XU  - upper bounds
2363b242c63SJacob Faibussowitsch . T1  - work vector
2373b242c63SJacob Faibussowitsch - T2  - work vector
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay   Output Parameters:
240a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result
241a7e14dcfSSatish Balay - Db - row scaling component of the result
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay   Level: developer
244a7e14dcfSSatish Balay 
245a1cb98faSBarry Smith .seealso: `Mat`, `VecFischer()`, `VecSFischer()`, `MatDSFischer()`
246a7e14dcfSSatish Balay @*/
247d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db)
248d71ae5a4SJacob Faibussowitsch {
249a7e14dcfSSatish Balay   PetscInt           i, nn;
25046bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u, *t2;
25146bdf8c8SLisandro Dalcin   PetscScalar       *da, *db, *t1;
252a7e14dcfSSatish Balay   PetscReal          ai, bi, ci, di, ei;
253a7e14dcfSSatish Balay 
254a7e14dcfSSatish Balay   PetscFunctionBegin;
2559566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nn));
2569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Con, &f));
2589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(XL, &l));
2599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(XU, &u));
2609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Da, &da));
2619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Db, &db));
2629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(T1, &t1));
2639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(T2, &t2));
264a7e14dcfSSatish Balay 
265a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
26646bdf8c8SLisandro Dalcin     da[i] = 0.0;
26746bdf8c8SLisandro Dalcin     db[i] = 0.0;
26846bdf8c8SLisandro Dalcin     t1[i] = 0.0;
269a7e14dcfSSatish Balay 
27046bdf8c8SLisandro Dalcin     if (PetscAbsScalar(f[i]) <= PETSC_MACHINE_EPSILON) {
27146bdf8c8SLisandro Dalcin       if (PetscRealPart(l[i]) > PETSC_NINFINITY && PetscAbsScalar(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) {
27246bdf8c8SLisandro Dalcin         t1[i] = 1.0;
27346bdf8c8SLisandro Dalcin         da[i] = 1.0;
274a7e14dcfSSatish Balay       }
275a7e14dcfSSatish Balay 
27646bdf8c8SLisandro Dalcin       if (PetscRealPart(u[i]) < PETSC_INFINITY && PetscAbsScalar(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) {
27746bdf8c8SLisandro Dalcin         t1[i] = 1.0;
27846bdf8c8SLisandro Dalcin         db[i] = 1.0;
279a7e14dcfSSatish Balay       }
280a7e14dcfSSatish Balay     }
281a7e14dcfSSatish Balay   }
282a7e14dcfSSatish Balay 
2839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(T1, &t1));
2849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(T2, &t2));
2859566063dSJacob Faibussowitsch   PetscCall(MatMult(jac, T1, T2));
2869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(T2, &t2));
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
28946bdf8c8SLisandro Dalcin     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) {
29046bdf8c8SLisandro Dalcin       da[i] = 0.0;
29146bdf8c8SLisandro Dalcin       db[i] = -1.0;
29246bdf8c8SLisandro Dalcin     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {
29346bdf8c8SLisandro Dalcin       if (PetscRealPart(db[i]) >= 1) {
294658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, PetscRealPart(t2[i]));
295a7e14dcfSSatish Balay 
29646bdf8c8SLisandro Dalcin         da[i] = -1.0 / ai - 1.0;
29746bdf8c8SLisandro Dalcin         db[i] = -t2[i] / ai - 1.0;
2982d0e5244SBarry Smith       } else {
299658c1fc4SLisandro Dalcin         bi = PetscRealPart(u[i]) - PetscRealPart(x[i]);
300658c1fc4SLisandro Dalcin         ai = fischnorm(bi, PetscRealPart(f[i]));
301a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
302a7e14dcfSSatish Balay 
30346bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
30446bdf8c8SLisandro Dalcin         db[i] = -f[i] / ai - 1.0;
305a7e14dcfSSatish Balay       }
30646bdf8c8SLisandro Dalcin     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) {
30746bdf8c8SLisandro Dalcin       if (PetscRealPart(da[i]) >= 1) {
308658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, PetscRealPart(t2[i]));
309a7e14dcfSSatish Balay 
31046bdf8c8SLisandro Dalcin         da[i] = 1.0 / ai - 1.0;
31146bdf8c8SLisandro Dalcin         db[i] = t2[i] / ai - 1.0;
3122d0e5244SBarry Smith       } else {
313658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(l[i]);
314658c1fc4SLisandro Dalcin         ai = fischnorm(bi, PetscRealPart(f[i]));
315a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
316a7e14dcfSSatish Balay 
31746bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
31846bdf8c8SLisandro Dalcin         db[i] = f[i] / ai - 1.0;
319a7e14dcfSSatish Balay       }
320658c1fc4SLisandro Dalcin     } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) {
32146bdf8c8SLisandro Dalcin       da[i] = -1.0;
32246bdf8c8SLisandro Dalcin       db[i] = 0.0;
3232d0e5244SBarry Smith     } else {
32446bdf8c8SLisandro Dalcin       if (PetscRealPart(db[i]) >= 1) {
325658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, PetscRealPart(t2[i]));
326a7e14dcfSSatish Balay 
32746bdf8c8SLisandro Dalcin         ci = 1.0 / ai + 1.0;
328658c1fc4SLisandro Dalcin         di = PetscRealPart(t2[i]) / ai + 1.0;
3292d0e5244SBarry Smith       } else {
330658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(u[i]);
331658c1fc4SLisandro Dalcin         ai = fischnorm(bi, PetscRealPart(f[i]));
332a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
333a7e14dcfSSatish Balay 
33446bdf8c8SLisandro Dalcin         ci = bi / ai + 1.0;
335658c1fc4SLisandro Dalcin         di = PetscRealPart(f[i]) / ai + 1.0;
336a7e14dcfSSatish Balay       }
337a7e14dcfSSatish Balay 
33846bdf8c8SLisandro Dalcin       if (PetscRealPart(da[i]) >= 1) {
339658c1fc4SLisandro Dalcin         bi = ci + di * PetscRealPart(t2[i]);
340658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, bi);
341a7e14dcfSSatish Balay 
34246bdf8c8SLisandro Dalcin         bi = bi / ai - 1.0;
34346bdf8c8SLisandro Dalcin         ai = 1.0 / ai - 1.0;
3442d0e5244SBarry Smith       } else {
345658c1fc4SLisandro Dalcin         ei = Fischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]));
346658c1fc4SLisandro Dalcin         ai = fischnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei);
347a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
348a7e14dcfSSatish Balay 
34946bdf8c8SLisandro Dalcin         bi = ei / ai - 1.0;
350658c1fc4SLisandro Dalcin         ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0;
351a7e14dcfSSatish Balay       }
352a7e14dcfSSatish Balay 
353a7e14dcfSSatish Balay       da[i] = ai + bi * ci;
354a7e14dcfSSatish Balay       db[i] = bi * di;
355a7e14dcfSSatish Balay     }
356a7e14dcfSSatish Balay   }
357a7e14dcfSSatish Balay 
3589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Da, &da));
3599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Db, &db));
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Con, &f));
3629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(XL, &l));
3639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(XU, &u));
3649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(T2, &t2));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3668e3154b5SSatish Balay }
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay /*@
369235fd6e6SBarry Smith   MatDSFischer - Calculates an element of the B-subdifferential of the
370a7e14dcfSSatish Balay   smoothed Fischer-Burmeister function for complementarity problems.
371a7e14dcfSSatish Balay 
372c3339decSBarry Smith   Collective
373a7e14dcfSSatish Balay 
374a7e14dcfSSatish Balay   Input Parameters:
375a7e14dcfSSatish Balay + jac - the jacobian of f at X
376a7e14dcfSSatish Balay . X   - current point
377e056e8ceSJacob Faibussowitsch . Con - constraint function evaluated at X
378a7e14dcfSSatish Balay . XL  - lower bounds
379a7e14dcfSSatish Balay . XU  - upper bounds
380a7e14dcfSSatish Balay . mu  - smoothing parameter
381a7e14dcfSSatish Balay . T1  - work vector
382a7e14dcfSSatish Balay - T2  - work vector
383a7e14dcfSSatish Balay 
384d8d19677SJose E. Roman   Output Parameters:
385a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result
386a7e14dcfSSatish Balay . Db - row scaling component of the result
387a7e14dcfSSatish Balay - Dm - derivative with respect to scaling parameter
388a7e14dcfSSatish Balay 
389a7e14dcfSSatish Balay   Level: developer
390a7e14dcfSSatish Balay 
391a1cb98faSBarry Smith .seealso: `Mat`, `VecFischer()`, `VecDFischer()`, `MatDFischer()`
392a7e14dcfSSatish Balay @*/
393d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDSFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, PetscReal mu, Vec T1, Vec T2, Vec Da, Vec Db, Vec Dm)
394d71ae5a4SJacob Faibussowitsch {
395a7e14dcfSSatish Balay   PetscInt           i, nn;
39646bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
39746bdf8c8SLisandro Dalcin   PetscScalar       *da, *db, *dm;
398a7e14dcfSSatish Balay   PetscReal          ai, bi, ci, di, ei, fi;
399a7e14dcfSSatish Balay 
400a7e14dcfSSatish Balay   PetscFunctionBegin;
401a7e14dcfSSatish Balay   if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) {
4029566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Dm));
4039566063dSJacob Faibussowitsch     PetscCall(MatDFischer(jac, X, Con, XL, XU, T1, T2, Da, Db));
4042d0e5244SBarry Smith   } else {
4059566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(X, &nn));
4069566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
4079566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Con, &f));
4089566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &l));
4099566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &u));
4109566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Da, &da));
4119566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Db, &db));
4129566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Dm, &dm));
413a7e14dcfSSatish Balay 
414a7e14dcfSSatish Balay     for (i = 0; i < nn; ++i) {
41546bdf8c8SLisandro Dalcin       if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) {
416a7e14dcfSSatish Balay         da[i] = -mu;
41746bdf8c8SLisandro Dalcin         db[i] = -1.0;
418a7e14dcfSSatish Balay         dm[i] = -x[i];
41946bdf8c8SLisandro Dalcin       } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {
420658c1fc4SLisandro Dalcin         bi = PetscRealPart(u[i]) - PetscRealPart(x[i]);
421658c1fc4SLisandro Dalcin         ai = fischsnorm(bi, PetscRealPart(f[i]), mu);
422a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
423a7e14dcfSSatish Balay 
42446bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
425658c1fc4SLisandro Dalcin         db[i] = -PetscRealPart(f[i]) / ai - 1.0;
426a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
42746bdf8c8SLisandro Dalcin       } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) {
428658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(l[i]);
429658c1fc4SLisandro Dalcin         ai = fischsnorm(bi, PetscRealPart(f[i]), mu);
430a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
431a7e14dcfSSatish Balay 
43246bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
433658c1fc4SLisandro Dalcin         db[i] = PetscRealPart(f[i]) / ai - 1.0;
434a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
435658c1fc4SLisandro Dalcin       } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) {
43646bdf8c8SLisandro Dalcin         da[i] = -1.0;
43746bdf8c8SLisandro Dalcin         db[i] = 0.0;
43846bdf8c8SLisandro Dalcin         dm[i] = 0.0;
4392d0e5244SBarry Smith       } else {
440658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(u[i]);
441658c1fc4SLisandro Dalcin         ai = fischsnorm(bi, PetscRealPart(f[i]), mu);
442a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
443a7e14dcfSSatish Balay 
44446bdf8c8SLisandro Dalcin         ci = bi / ai + 1.0;
445658c1fc4SLisandro Dalcin         di = PetscRealPart(f[i]) / ai + 1.0;
446a7e14dcfSSatish Balay         fi = 2.0 * mu / ai;
447a7e14dcfSSatish Balay 
448658c1fc4SLisandro Dalcin         ei = SFischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]), mu);
449658c1fc4SLisandro Dalcin         ai = fischsnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei, mu);
450a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
451a7e14dcfSSatish Balay 
45246bdf8c8SLisandro Dalcin         bi = ei / ai - 1.0;
453a7e14dcfSSatish Balay         ei = 2.0 * mu / ei;
454658c1fc4SLisandro Dalcin         ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0;
455a7e14dcfSSatish Balay 
456a7e14dcfSSatish Balay         da[i] = ai + bi * ci;
457a7e14dcfSSatish Balay         db[i] = bi * di;
458a7e14dcfSSatish Balay         dm[i] = ei + bi * fi;
459a7e14dcfSSatish Balay       }
460a7e14dcfSSatish Balay     }
461a7e14dcfSSatish Balay 
4629566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
4639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Con, &f));
4649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &l));
4659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &u));
4669566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Da, &da));
4679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Db, &db));
4689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Dm, &dm));
469a7e14dcfSSatish Balay   }
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471a7e14dcfSSatish Balay }
472a7e14dcfSSatish Balay 
473d71ae5a4SJacob Faibussowitsch static inline PetscReal ST_InternalPN(PetscScalar in, PetscReal lb, PetscReal ub)
474d71ae5a4SJacob Faibussowitsch {
475835f2295SStefano Zampini   return PetscMax(0, PetscRealPart(in) - ub) - PetscMax(0, -PetscRealPart(in) - PetscAbsReal(lb));
4768370d7cdSHansol Suh }
4778370d7cdSHansol Suh 
478d71ae5a4SJacob Faibussowitsch static inline PetscReal ST_InternalNN(PetscScalar in, PetscReal lb, PetscReal ub)
479d71ae5a4SJacob Faibussowitsch {
480835f2295SStefano Zampini   return PetscMax(0, PetscRealPart(in) + PetscAbsReal(ub)) - PetscMax(0, -PetscRealPart(in) - PetscAbsReal(lb));
4818370d7cdSHansol Suh }
4828370d7cdSHansol Suh 
483d71ae5a4SJacob Faibussowitsch static inline PetscReal ST_InternalPP(PetscScalar in, PetscReal lb, PetscReal ub)
484d71ae5a4SJacob Faibussowitsch {
485835f2295SStefano Zampini   return PetscMax(0, PetscRealPart(in) - ub) + PetscMin(0, PetscRealPart(in) - lb);
4868370d7cdSHansol Suh }
4878370d7cdSHansol Suh 
4888370d7cdSHansol Suh /*@
4898370d7cdSHansol Suh   TaoSoftThreshold - Calculates soft thresholding routine with input vector
4908370d7cdSHansol Suh   and given lower and upper bound and returns it to output vector.
4918370d7cdSHansol Suh 
4928370d7cdSHansol Suh   Input Parameters:
4938370d7cdSHansol Suh + in - input vector to be thresholded
4948370d7cdSHansol Suh . lb - lower bound
495f0fc11ceSJed Brown - ub - upper bound
4968370d7cdSHansol Suh 
49797bb3fdcSJose E. Roman   Output Parameter:
4988370d7cdSHansol Suh . out - Soft thresholded output vector
4998370d7cdSHansol Suh 
5008370d7cdSHansol Suh   Notes:
5018370d7cdSHansol Suh   Soft thresholding is defined as
5028370d7cdSHansol Suh   \[ S(input,lb,ub) =
5038370d7cdSHansol Suh   \begin{cases}
5048370d7cdSHansol Suh   input - ub  \text{input > ub} \\
5058370d7cdSHansol Suh   0           \text{lb =< input <= ub} \\
5068370d7cdSHansol Suh   input + lb  \text{input < lb} \\
5078370d7cdSHansol Suh   \]
5088370d7cdSHansol Suh 
5098370d7cdSHansol Suh   Level: developer
5108370d7cdSHansol Suh 
511a1cb98faSBarry Smith .seealso: `Tao`, `Vec`
5128370d7cdSHansol Suh @*/
513d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSoftThreshold(Vec in, PetscReal lb, PetscReal ub, Vec out)
514d71ae5a4SJacob Faibussowitsch {
5158370d7cdSHansol Suh   PetscInt     i, nlocal, mlocal;
5168370d7cdSHansol Suh   PetscScalar *inarray, *outarray;
5178370d7cdSHansol Suh 
5188370d7cdSHansol Suh   PetscFunctionBegin;
5199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayPair(in, out, &inarray, &outarray));
5209566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(in, &nlocal));
52184430a0dSHansol Suh   PetscCall(VecGetLocalSize(out, &mlocal));
5228370d7cdSHansol Suh 
5233c859ba3SBarry Smith   PetscCheck(nlocal == mlocal, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input and output vectors need to be of same size.");
52484430a0dSHansol Suh   if (lb == ub) {
52584430a0dSHansol Suh     PetscCall(VecRestoreArrayPair(in, out, &inarray, &outarray));
52684430a0dSHansol Suh     PetscFunctionReturn(PETSC_SUCCESS);
52784430a0dSHansol Suh   }
5283c859ba3SBarry Smith   PetscCheck(lb < ub, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Lower bound needs to be lower than upper bound.");
5298370d7cdSHansol Suh 
5308370d7cdSHansol Suh   if (ub >= 0 && lb < 0) {
5318370d7cdSHansol Suh     for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalPN(inarray[i], lb, ub);
5328370d7cdSHansol Suh   } else if (ub < 0 && lb < 0) {
5338370d7cdSHansol Suh     for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalNN(inarray[i], lb, ub);
5348370d7cdSHansol Suh   } else {
5358370d7cdSHansol Suh     for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalPP(inarray[i], lb, ub);
5368370d7cdSHansol Suh   }
5378370d7cdSHansol Suh 
5389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayPair(in, out, &inarray, &outarray));
5393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5408370d7cdSHansol Suh }
541