xref: /petsc/src/tao/util/tao_util.c (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2ba92ff59SBarry Smith #include <petsctao.h>      /*I "petsctao.h" I*/
38370d7cdSHansol Suh #include <petscsys.h>
4a7e14dcfSSatish Balay 
59fbee547SJacob Faibussowitsch static inline PetscReal Fischer(PetscReal a, PetscReal b)
6a7e14dcfSSatish Balay {
7a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
8a7e14dcfSSatish Balay    if (a + b <= 0) {
946bdf8c8SLisandro Dalcin      return PetscSqrtReal(a*a + b*b) - (a + b);
10a7e14dcfSSatish Balay    }
1146bdf8c8SLisandro Dalcin    return -2.0*a*b / (PetscSqrtReal(a*a + b*b) + (a + b));
12a7e14dcfSSatish Balay }
13a7e14dcfSSatish Balay 
14a7e14dcfSSatish Balay /*@
15a7e14dcfSSatish Balay    VecFischer - Evaluates the Fischer-Burmeister function for complementarity
16a7e14dcfSSatish Balay    problems.
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay    Logically Collective on vectors
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay    Input Parameters:
21a7e14dcfSSatish Balay +  X - current point
22a7e14dcfSSatish Balay .  F - function evaluated at x
23a7e14dcfSSatish Balay .  L - lower bounds
24a7e14dcfSSatish Balay -  U - upper bounds
25a7e14dcfSSatish Balay 
26f899ff85SJose E. Roman    Output Parameter:
27a7e14dcfSSatish Balay .  FB - The Fischer-Burmeister function vector
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay    Notes:
30a7e14dcfSSatish Balay    The Fischer-Burmeister function is defined as
31a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b) - a - b
32a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
33a7e14dcfSSatish Balay    system of equations.
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay    The result of this function is done by cases:
36a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i]
37a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i])
38a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i])
39a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u]))
40a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay    Level: developer
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay @*/
45a7e14dcfSSatish Balay PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB)
46a7e14dcfSSatish Balay {
4746bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
4846bdf8c8SLisandro Dalcin   PetscScalar       *fb;
49a7e14dcfSSatish Balay   PetscReal         xval, fval, lval, uval;
50a7e14dcfSSatish Balay   PetscInt          low[5], high[5], n, i;
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay   PetscFunctionBegin;
53a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,1);
54a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID,2);
55*76be6f4fSStefano Zampini   if (L) PetscValidHeaderSpecific(L, VEC_CLASSID,3);
56*76be6f4fSStefano Zampini   if (U) PetscValidHeaderSpecific(U, VEC_CLASSID,4);
57064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(FB, VEC_CLASSID,5);
58a7e14dcfSSatish Balay 
59*76be6f4fSStefano Zampini   if (!L && !U) {
60*76be6f4fSStefano Zampini     PetscCall(VecAXPBY(FB,-1.0,0.0,F));
61*76be6f4fSStefano Zampini     PetscFunctionReturn(0);
62*76be6f4fSStefano Zampini   }
63*76be6f4fSStefano Zampini 
649566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, low, high));
659566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(F, low + 1, high + 1));
669566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(L, low + 2, high + 2));
679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(U, low + 3, high + 3));
689566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4));
69a7e14dcfSSatish Balay 
70a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
713c859ba3SBarry Smith     PetscCheck(low[0] == low[i] && high[0] == high[i],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vectors must be identically loaded over processors");
72a7e14dcfSSatish Balay   }
73a7e14dcfSSatish Balay 
749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(L, &l));
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(FB, &fb));
79a7e14dcfSSatish Balay 
809566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
83*76be6f4fSStefano Zampini     xval = PetscRealPart(x[i]);
84*76be6f4fSStefano Zampini     fval = PetscRealPart(f[i]);
85*76be6f4fSStefano Zampini     lval = PetscRealPart(l[i]);
86*76be6f4fSStefano Zampini     uval = PetscRealPart(u[i]);
87a7e14dcfSSatish Balay 
88*76be6f4fSStefano Zampini     if (lval <= -PETSC_INFINITY && uval >= PETSC_INFINITY) {
89a7e14dcfSSatish Balay       fb[i] = -fval;
90e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
91a7e14dcfSSatish Balay       fb[i] = -Fischer(uval - xval, -fval);
92e270355aSBarry Smith     } else if (uval >=  PETSC_INFINITY) {
93a7e14dcfSSatish Balay       fb[i] = Fischer(xval - lval,  fval);
942d0e5244SBarry Smith     } else if (lval == uval) {
95a7e14dcfSSatish Balay       fb[i] = lval - xval;
962d0e5244SBarry Smith     } else {
97a7e14dcfSSatish Balay       fval  = Fischer(uval - xval, -fval);
98a7e14dcfSSatish Balay       fb[i] = Fischer(xval - lval,  fval);
99a7e14dcfSSatish Balay     }
100a7e14dcfSSatish Balay   }
101a7e14dcfSSatish Balay 
1029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
1049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(L, &l));
1059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(FB, &fb));
107a7e14dcfSSatish Balay   PetscFunctionReturn(0);
108a7e14dcfSSatish Balay }
109a7e14dcfSSatish Balay 
1109fbee547SJacob Faibussowitsch static inline PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c)
111a7e14dcfSSatish Balay {
112a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
113a7e14dcfSSatish Balay    if (a + b <= 0) {
1143f6ba705SLisandro Dalcin      return PetscSqrtReal(a*a + b*b + 2.0*c*c) - (a + b);
115a7e14dcfSSatish Balay    }
1163f6ba705SLisandro Dalcin    return 2.0*(c*c - a*b) / (PetscSqrtReal(a*a + b*b + 2.0*c*c) + (a + b));
117a7e14dcfSSatish Balay }
118a7e14dcfSSatish Balay 
119a7e14dcfSSatish Balay /*@
120a7e14dcfSSatish Balay    VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for
121a7e14dcfSSatish Balay    complementarity problems.
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay    Logically Collective on vectors
124a7e14dcfSSatish Balay 
125a7e14dcfSSatish Balay    Input Parameters:
126a7e14dcfSSatish Balay +  X - current point
127a7e14dcfSSatish Balay .  F - function evaluated at x
128a7e14dcfSSatish Balay .  L - lower bounds
129a7e14dcfSSatish Balay .  U - upper bounds
130a7e14dcfSSatish Balay -  mu - smoothing parameter
131a7e14dcfSSatish Balay 
132f899ff85SJose E. Roman    Output Parameter:
133a7e14dcfSSatish Balay .  FB - The Smoothed Fischer-Burmeister function vector
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay    Notes:
136a7e14dcfSSatish Balay    The Smoothed Fischer-Burmeister function is defined as
137a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b
138a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
139a7e14dcfSSatish Balay    system of equations.
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay    The result of this function is done by cases:
142a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i] - 2*mu*x[i]
143a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i], mu)
144a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i], mu)
145a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu)
146a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay    Level: developer
149a7e14dcfSSatish Balay 
150db781477SPatrick Sanan .seealso `VecFischer()`
151a7e14dcfSSatish Balay @*/
152a7e14dcfSSatish Balay PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB)
153a7e14dcfSSatish Balay {
15446bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
15546bdf8c8SLisandro Dalcin   PetscScalar       *fb;
156a7e14dcfSSatish Balay   PetscReal         xval, fval, lval, uval;
157a7e14dcfSSatish Balay   PetscInt          low[5], high[5], n, i;
158a7e14dcfSSatish Balay 
159a7e14dcfSSatish Balay   PetscFunctionBegin;
160a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,1);
161a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID,2);
162a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID,3);
163a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID,4);
164a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID,6);
165a7e14dcfSSatish Balay 
1669566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, low, high));
1679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(F, low + 1, high + 1));
1689566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(L, low + 2, high + 2));
1699566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(U, low + 3, high + 3));
1709566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4));
171a7e14dcfSSatish Balay 
172a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
1733c859ba3SBarry Smith     PetscCheck(low[0] == low[i] && high[0] == high[i],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vectors must be identically loaded over processors");
174a7e14dcfSSatish Balay   }
175a7e14dcfSSatish Balay 
1769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(L, &l));
1799566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(FB, &fb));
181a7e14dcfSSatish Balay 
1829566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
183a7e14dcfSSatish Balay 
184a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
185658c1fc4SLisandro Dalcin     xval = PetscRealPart(*x++); fval = PetscRealPart(*f++);
186658c1fc4SLisandro Dalcin     lval = PetscRealPart(*l++); uval = PetscRealPart(*u++);
187a7e14dcfSSatish Balay 
188e270355aSBarry Smith     if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) {
189a7e14dcfSSatish Balay       (*fb++) = -fval - mu*xval;
190e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
191a7e14dcfSSatish Balay       (*fb++) = -SFischer(uval - xval, -fval, mu);
192e270355aSBarry Smith     } else if (uval >=  PETSC_INFINITY) {
193a7e14dcfSSatish Balay       (*fb++) =  SFischer(xval - lval,  fval, mu);
1942d0e5244SBarry Smith     } else if (lval == uval) {
195a7e14dcfSSatish Balay       (*fb++) = lval - xval;
1962d0e5244SBarry Smith     } else {
197a7e14dcfSSatish Balay       fval    =  SFischer(uval - xval, -fval, mu);
198a7e14dcfSSatish Balay       (*fb++) =  SFischer(xval - lval,  fval, mu);
199a7e14dcfSSatish Balay     }
200a7e14dcfSSatish Balay   }
201a7e14dcfSSatish Balay   x -= n; f -= n; l -=n; u -= n; fb -= n;
202a7e14dcfSSatish Balay 
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
2059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(L, &l));
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(FB, &fb));
208a7e14dcfSSatish Balay   PetscFunctionReturn(0);
209a7e14dcfSSatish Balay }
210a7e14dcfSSatish Balay 
2119fbee547SJacob Faibussowitsch static inline PetscReal fischnorm(PetscReal a, PetscReal b)
212a7e14dcfSSatish Balay {
213658c1fc4SLisandro Dalcin   return PetscSqrtReal(a*a + b*b);
214a7e14dcfSSatish Balay }
215a7e14dcfSSatish Balay 
2169fbee547SJacob Faibussowitsch static inline PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c)
217a7e14dcfSSatish Balay {
218658c1fc4SLisandro Dalcin   return PetscSqrtReal(a*a + b*b + 2.0*c*c);
219a7e14dcfSSatish Balay }
220a7e14dcfSSatish Balay 
221a7e14dcfSSatish Balay /*@
222235fd6e6SBarry Smith    MatDFischer - Calculates an element of the B-subdifferential of the
223a7e14dcfSSatish Balay    Fischer-Burmeister function for complementarity problems.
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay    Collective on jac
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay    Input Parameters:
228a7e14dcfSSatish Balay +  jac - the jacobian of f at X
229a7e14dcfSSatish Balay .  X - current point
230a7e14dcfSSatish Balay .  Con - constraints function evaluated at X
231a7e14dcfSSatish Balay .  XL - lower bounds
232a7e14dcfSSatish Balay .  XU - upper bounds
233a7e14dcfSSatish Balay .  t1 - work vector
234a7e14dcfSSatish Balay -  t2 - work vector
235a7e14dcfSSatish Balay 
236a7e14dcfSSatish Balay    Output Parameters:
237a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
238a7e14dcfSSatish Balay -  Db - row scaling component of the result
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay    Level: developer
241a7e14dcfSSatish Balay 
242db781477SPatrick Sanan .seealso: `VecFischer()`
243a7e14dcfSSatish Balay @*/
244235fd6e6SBarry Smith PetscErrorCode MatDFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db)
245a7e14dcfSSatish Balay {
246a7e14dcfSSatish Balay   PetscInt          i,nn;
24746bdf8c8SLisandro Dalcin   const PetscScalar *x,*f,*l,*u,*t2;
24846bdf8c8SLisandro Dalcin   PetscScalar       *da,*db,*t1;
249a7e14dcfSSatish Balay   PetscReal          ai,bi,ci,di,ei;
250a7e14dcfSSatish Balay 
251a7e14dcfSSatish Balay   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&nn));
2539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Con,&f));
2559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(XL,&l));
2569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(XU,&u));
2579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Da,&da));
2589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Db,&db));
2599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(T1,&t1));
2609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(T2,&t2));
261a7e14dcfSSatish Balay 
262a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
26346bdf8c8SLisandro Dalcin     da[i] = 0.0;
26446bdf8c8SLisandro Dalcin     db[i] = 0.0;
26546bdf8c8SLisandro Dalcin     t1[i] = 0.0;
266a7e14dcfSSatish Balay 
26746bdf8c8SLisandro Dalcin     if (PetscAbsScalar(f[i]) <= PETSC_MACHINE_EPSILON) {
26846bdf8c8SLisandro Dalcin       if (PetscRealPart(l[i]) > PETSC_NINFINITY && PetscAbsScalar(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) {
26946bdf8c8SLisandro Dalcin         t1[i] = 1.0;
27046bdf8c8SLisandro Dalcin         da[i] = 1.0;
271a7e14dcfSSatish Balay       }
272a7e14dcfSSatish Balay 
27346bdf8c8SLisandro Dalcin       if (PetscRealPart(u[i]) <  PETSC_INFINITY && PetscAbsScalar(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) {
27446bdf8c8SLisandro Dalcin         t1[i] = 1.0;
27546bdf8c8SLisandro Dalcin         db[i] = 1.0;
276a7e14dcfSSatish Balay       }
277a7e14dcfSSatish Balay     }
278a7e14dcfSSatish Balay   }
279a7e14dcfSSatish Balay 
2809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(T1,&t1));
2819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(T2,&t2));
2829566063dSJacob Faibussowitsch   PetscCall(MatMult(jac,T1,T2));
2839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(T2,&t2));
284a7e14dcfSSatish Balay 
285a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
28646bdf8c8SLisandro Dalcin     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) {
28746bdf8c8SLisandro Dalcin       da[i] = 0.0;
28846bdf8c8SLisandro Dalcin       db[i] = -1.0;
28946bdf8c8SLisandro Dalcin     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {
29046bdf8c8SLisandro Dalcin       if (PetscRealPart(db[i]) >= 1) {
291658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, PetscRealPart(t2[i]));
292a7e14dcfSSatish Balay 
29346bdf8c8SLisandro Dalcin         da[i] = -1.0 / ai - 1.0;
29446bdf8c8SLisandro Dalcin         db[i] = -t2[i] / ai - 1.0;
2952d0e5244SBarry Smith       } else {
296658c1fc4SLisandro Dalcin         bi = PetscRealPart(u[i]) - PetscRealPart(x[i]);
297658c1fc4SLisandro Dalcin         ai = fischnorm(bi, PetscRealPart(f[i]));
298a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
299a7e14dcfSSatish Balay 
30046bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
30146bdf8c8SLisandro Dalcin         db[i] = -f[i] / ai - 1.0;
302a7e14dcfSSatish Balay       }
30346bdf8c8SLisandro Dalcin     } else if (PetscRealPart(u[i]) >=  PETSC_INFINITY) {
30446bdf8c8SLisandro Dalcin       if (PetscRealPart(da[i]) >= 1) {
305658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, PetscRealPart(t2[i]));
306a7e14dcfSSatish Balay 
30746bdf8c8SLisandro Dalcin         da[i] = 1.0 / ai - 1.0;
30846bdf8c8SLisandro Dalcin         db[i] = t2[i] / ai - 1.0;
3092d0e5244SBarry Smith       } else {
310658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(l[i]);
311658c1fc4SLisandro Dalcin         ai = fischnorm(bi, PetscRealPart(f[i]));
312a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
313a7e14dcfSSatish Balay 
31446bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
31546bdf8c8SLisandro Dalcin         db[i] = f[i] / ai - 1.0;
316a7e14dcfSSatish Balay       }
317658c1fc4SLisandro Dalcin     } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) {
31846bdf8c8SLisandro Dalcin       da[i] = -1.0;
31946bdf8c8SLisandro Dalcin       db[i] = 0.0;
3202d0e5244SBarry Smith     } else {
32146bdf8c8SLisandro Dalcin       if (PetscRealPart(db[i]) >= 1) {
322658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, PetscRealPart(t2[i]));
323a7e14dcfSSatish Balay 
32446bdf8c8SLisandro Dalcin         ci = 1.0 / ai + 1.0;
325658c1fc4SLisandro Dalcin         di = PetscRealPart(t2[i]) / ai + 1.0;
3262d0e5244SBarry Smith       } else {
327658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(u[i]);
328658c1fc4SLisandro Dalcin         ai = fischnorm(bi, PetscRealPart(f[i]));
329a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
330a7e14dcfSSatish Balay 
33146bdf8c8SLisandro Dalcin         ci = bi / ai + 1.0;
332658c1fc4SLisandro Dalcin         di = PetscRealPart(f[i]) / ai + 1.0;
333a7e14dcfSSatish Balay       }
334a7e14dcfSSatish Balay 
33546bdf8c8SLisandro Dalcin       if (PetscRealPart(da[i]) >= 1) {
336658c1fc4SLisandro Dalcin         bi = ci + di*PetscRealPart(t2[i]);
337658c1fc4SLisandro Dalcin         ai = fischnorm(1.0, bi);
338a7e14dcfSSatish Balay 
33946bdf8c8SLisandro Dalcin         bi = bi / ai - 1.0;
34046bdf8c8SLisandro Dalcin         ai = 1.0 / ai - 1.0;
3412d0e5244SBarry Smith       } else {
342658c1fc4SLisandro Dalcin         ei = Fischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]));
343658c1fc4SLisandro Dalcin         ai = fischnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei);
344a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
345a7e14dcfSSatish Balay 
34646bdf8c8SLisandro Dalcin         bi = ei / ai - 1.0;
347658c1fc4SLisandro Dalcin         ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0;
348a7e14dcfSSatish Balay       }
349a7e14dcfSSatish Balay 
350a7e14dcfSSatish Balay       da[i] = ai + bi*ci;
351a7e14dcfSSatish Balay       db[i] = bi*di;
352a7e14dcfSSatish Balay     }
353a7e14dcfSSatish Balay   }
354a7e14dcfSSatish Balay 
3559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Da,&da));
3569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Db,&db));
3579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
3589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Con,&f));
3599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(XL,&l));
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(XU,&u));
3619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(T2,&t2));
362a7e14dcfSSatish Balay   PetscFunctionReturn(0);
3638e3154b5SSatish Balay }
364a7e14dcfSSatish Balay 
365a7e14dcfSSatish Balay /*@
366235fd6e6SBarry Smith    MatDSFischer - Calculates an element of the B-subdifferential of the
367a7e14dcfSSatish Balay    smoothed Fischer-Burmeister function for complementarity problems.
368a7e14dcfSSatish Balay 
369a7e14dcfSSatish Balay    Collective on jac
370a7e14dcfSSatish Balay 
371a7e14dcfSSatish Balay    Input Parameters:
372a7e14dcfSSatish Balay +  jac - the jacobian of f at X
373a7e14dcfSSatish Balay .  X - current point
374a7e14dcfSSatish Balay .  F - constraint function evaluated at X
375a7e14dcfSSatish Balay .  XL - lower bounds
376a7e14dcfSSatish Balay .  XU - upper bounds
377a7e14dcfSSatish Balay .  mu - smoothing parameter
378a7e14dcfSSatish Balay .  T1 - work vector
379a7e14dcfSSatish Balay -  T2 - work vector
380a7e14dcfSSatish Balay 
381d8d19677SJose E. Roman    Output Parameters:
382a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
383a7e14dcfSSatish Balay .  Db - row scaling component of the result
384a7e14dcfSSatish Balay -  Dm - derivative with respect to scaling parameter
385a7e14dcfSSatish Balay 
386a7e14dcfSSatish Balay    Level: developer
387a7e14dcfSSatish Balay 
388db781477SPatrick Sanan .seealso `MatDFischer()`
389a7e14dcfSSatish Balay @*/
390235fd6e6SBarry Smith PetscErrorCode MatDSFischer(Mat jac, Vec X, Vec Con,Vec XL, Vec XU, PetscReal mu,Vec T1, Vec T2,Vec Da, Vec Db, Vec Dm)
391a7e14dcfSSatish Balay {
392a7e14dcfSSatish Balay   PetscInt          i,nn;
39346bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
39446bdf8c8SLisandro Dalcin   PetscScalar       *da, *db, *dm;
395a7e14dcfSSatish Balay   PetscReal         ai, bi, ci, di, ei, fi;
396a7e14dcfSSatish Balay 
397a7e14dcfSSatish Balay   PetscFunctionBegin;
398a7e14dcfSSatish Balay   if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) {
3999566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Dm));
4009566063dSJacob Faibussowitsch     PetscCall(MatDFischer(jac, X, Con, XL, XU, T1, T2, Da, Db));
4012d0e5244SBarry Smith   } else {
4029566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(X,&nn));
4039566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X,&x));
4049566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Con,&f));
4059566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL,&l));
4069566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU,&u));
4079566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Da,&da));
4089566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Db,&db));
4099566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Dm,&dm));
410a7e14dcfSSatish Balay 
411a7e14dcfSSatish Balay     for (i = 0; i < nn; ++i) {
41246bdf8c8SLisandro Dalcin       if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) {
413a7e14dcfSSatish Balay         da[i] = -mu;
41446bdf8c8SLisandro Dalcin         db[i] = -1.0;
415a7e14dcfSSatish Balay         dm[i] = -x[i];
41646bdf8c8SLisandro Dalcin       } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {
417658c1fc4SLisandro Dalcin         bi = PetscRealPart(u[i]) - PetscRealPart(x[i]);
418658c1fc4SLisandro Dalcin         ai = fischsnorm(bi, PetscRealPart(f[i]), mu);
419a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
420a7e14dcfSSatish Balay 
42146bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
422658c1fc4SLisandro Dalcin         db[i] = -PetscRealPart(f[i]) / ai - 1.0;
423a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
42446bdf8c8SLisandro Dalcin       } else if (PetscRealPart(u[i]) >=  PETSC_INFINITY) {
425658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(l[i]);
426658c1fc4SLisandro Dalcin         ai = fischsnorm(bi, PetscRealPart(f[i]), mu);
427a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
428a7e14dcfSSatish Balay 
42946bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
430658c1fc4SLisandro Dalcin         db[i] = PetscRealPart(f[i]) / ai - 1.0;
431a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
432658c1fc4SLisandro Dalcin       } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) {
43346bdf8c8SLisandro Dalcin         da[i] = -1.0;
43446bdf8c8SLisandro Dalcin         db[i] = 0.0;
43546bdf8c8SLisandro Dalcin         dm[i] = 0.0;
4362d0e5244SBarry Smith       } else {
437658c1fc4SLisandro Dalcin         bi = PetscRealPart(x[i]) - PetscRealPart(u[i]);
438658c1fc4SLisandro Dalcin         ai = fischsnorm(bi, PetscRealPart(f[i]), mu);
439a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
440a7e14dcfSSatish Balay 
44146bdf8c8SLisandro Dalcin         ci = bi / ai + 1.0;
442658c1fc4SLisandro Dalcin         di = PetscRealPart(f[i]) / ai + 1.0;
443a7e14dcfSSatish Balay         fi = 2.0 * mu / ai;
444a7e14dcfSSatish Balay 
445658c1fc4SLisandro Dalcin         ei = SFischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]), mu);
446658c1fc4SLisandro Dalcin         ai = fischsnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei, mu);
447a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
448a7e14dcfSSatish Balay 
44946bdf8c8SLisandro Dalcin         bi = ei / ai - 1.0;
450a7e14dcfSSatish Balay         ei = 2.0 * mu / ei;
451658c1fc4SLisandro Dalcin         ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0;
452a7e14dcfSSatish Balay 
453a7e14dcfSSatish Balay         da[i] = ai + bi*ci;
454a7e14dcfSSatish Balay         db[i] = bi*di;
455a7e14dcfSSatish Balay         dm[i] = ei + bi*fi;
456a7e14dcfSSatish Balay       }
457a7e14dcfSSatish Balay     }
458a7e14dcfSSatish Balay 
4599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X,&x));
4609566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Con,&f));
4619566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL,&l));
4629566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU,&u));
4639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Da,&da));
4649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Db,&db));
4659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Dm,&dm));
466a7e14dcfSSatish Balay   }
467a7e14dcfSSatish Balay   PetscFunctionReturn(0);
468a7e14dcfSSatish Balay }
469a7e14dcfSSatish Balay 
4709fbee547SJacob Faibussowitsch static inline PetscReal ST_InternalPN(PetscScalar in, PetscReal lb, PetscReal ub)
4718370d7cdSHansol Suh {
4728370d7cdSHansol Suh   return PetscMax(0,(PetscReal)PetscRealPart(in)-ub) - PetscMax(0,-(PetscReal)PetscRealPart(in)-PetscAbsReal(lb));
4738370d7cdSHansol Suh }
4748370d7cdSHansol Suh 
4759fbee547SJacob Faibussowitsch static inline PetscReal ST_InternalNN(PetscScalar in, PetscReal lb, PetscReal ub)
4768370d7cdSHansol Suh {
4778370d7cdSHansol Suh   return PetscMax(0,(PetscReal)PetscRealPart(in) + PetscAbsReal(ub)) - PetscMax(0,-(PetscReal)PetscRealPart(in) - PetscAbsReal(lb));
4788370d7cdSHansol Suh }
4798370d7cdSHansol Suh 
4809fbee547SJacob Faibussowitsch static inline PetscReal ST_InternalPP(PetscScalar in, PetscReal lb, PetscReal ub)
4818370d7cdSHansol Suh {
4828370d7cdSHansol Suh   return PetscMax(0, (PetscReal)PetscRealPart(in)-ub) + PetscMin(0, (PetscReal)PetscRealPart(in) - lb);
4838370d7cdSHansol Suh }
4848370d7cdSHansol Suh 
4858370d7cdSHansol Suh /*@
4868370d7cdSHansol Suh    TaoSoftThreshold - Calculates soft thresholding routine with input vector
4878370d7cdSHansol Suh    and given lower and upper bound and returns it to output vector.
4888370d7cdSHansol Suh 
4898370d7cdSHansol Suh    Input Parameters:
4908370d7cdSHansol Suh +  in - input vector to be thresholded
4918370d7cdSHansol Suh .  lb - lower bound
492f0fc11ceSJed Brown -  ub - upper bound
4938370d7cdSHansol Suh 
49497bb3fdcSJose E. Roman    Output Parameter:
4958370d7cdSHansol Suh .  out - Soft thresholded output vector
4968370d7cdSHansol Suh 
4978370d7cdSHansol Suh    Notes:
4988370d7cdSHansol Suh    Soft thresholding is defined as
4998370d7cdSHansol Suh    \[ S(input,lb,ub) =
5008370d7cdSHansol Suh      \begin{cases}
5018370d7cdSHansol Suh     input - ub  \text{input > ub} \\
5028370d7cdSHansol Suh     0           \text{lb =< input <= ub} \\
5038370d7cdSHansol Suh     input + lb  \text{input < lb} \\
5048370d7cdSHansol Suh    \]
5058370d7cdSHansol Suh 
5068370d7cdSHansol Suh    Level: developer
5078370d7cdSHansol Suh 
5088370d7cdSHansol Suh @*/
5098370d7cdSHansol Suh PetscErrorCode TaoSoftThreshold(Vec in, PetscReal lb, PetscReal ub, Vec out)
5108370d7cdSHansol Suh {
5118370d7cdSHansol Suh   PetscInt       i, nlocal, mlocal;
5128370d7cdSHansol Suh   PetscScalar   *inarray, *outarray;
5138370d7cdSHansol Suh 
5148370d7cdSHansol Suh   PetscFunctionBegin;
5159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayPair(in, out, &inarray, &outarray));
5169566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(in, &nlocal));
5179566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(in, &mlocal));
5188370d7cdSHansol Suh 
5193c859ba3SBarry Smith   PetscCheck(nlocal == mlocal,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input and output vectors need to be of same size.");
5203c859ba3SBarry Smith   PetscCheck(lb < ub,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Lower bound needs to be lower than upper bound.");
5218370d7cdSHansol Suh 
5228370d7cdSHansol Suh   if (ub >= 0 && lb < 0) {
5238370d7cdSHansol Suh     for (i=0; i<nlocal; i++) outarray[i] = ST_InternalPN(inarray[i], lb, ub);
5248370d7cdSHansol Suh   } else if (ub < 0 && lb < 0) {
5258370d7cdSHansol Suh     for (i=0; i<nlocal; i++) outarray[i] = ST_InternalNN(inarray[i], lb, ub);
5268370d7cdSHansol Suh   } else {
5278370d7cdSHansol Suh     for (i=0; i<nlocal; i++) outarray[i] = ST_InternalPP(inarray[i], lb, ub);
5288370d7cdSHansol Suh   }
5298370d7cdSHansol Suh 
5309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayPair(in, out, &inarray, &outarray));
5318370d7cdSHansol Suh   PetscFunctionReturn(0);
5328370d7cdSHansol Suh }
533