xref: /petsc/src/tao/util/tao_util.c (revision ba92ff593176f3ffed64b48a0721b2817410e47a)
1aaa7dc30SBarry Smith #include <petsc-private/petscimpl.h>
2*ba92ff59SBarry Smith #include <petsctao.h>      /*I "petsctao.h" I*/
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal Fischer(PetscReal a, PetscReal b)
6a7e14dcfSSatish Balay {
7a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
8a7e14dcfSSatish Balay    if (a + b <= 0) {
9a7e14dcfSSatish Balay      return PetscSqrtScalar(a*a + b*b) - (a + b);
10a7e14dcfSSatish Balay    }
11a7e14dcfSSatish Balay    return -2.0*a*b / (PetscSqrtScalar(a*a + b*b) + (a + b));
12a7e14dcfSSatish Balay }
13a7e14dcfSSatish Balay 
14a7e14dcfSSatish Balay #undef __FUNCT__
15a7e14dcfSSatish Balay #define __FUNCT__ "VecFischer"
16a7e14dcfSSatish Balay /*@
17a7e14dcfSSatish Balay    VecFischer - Evaluates the Fischer-Burmeister function for complementarity
18a7e14dcfSSatish Balay    problems.
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay    Logically Collective on vectors
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay    Input Parameters:
23a7e14dcfSSatish Balay +  X - current point
24a7e14dcfSSatish Balay .  F - function evaluated at x
25a7e14dcfSSatish Balay .  L - lower bounds
26a7e14dcfSSatish Balay -  U - upper bounds
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay    Output Parameters:
29a7e14dcfSSatish Balay .  FB - The Fischer-Burmeister function vector
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay    Notes:
32a7e14dcfSSatish Balay    The Fischer-Burmeister function is defined as
33a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b) - a - b
34a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
35a7e14dcfSSatish Balay    system of equations.
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay    The result of this function is done by cases:
38a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i]
39a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i])
40a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i])
41a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u]))
42a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay    Level: developer
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay @*/
47a7e14dcfSSatish Balay PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB)
48a7e14dcfSSatish Balay {
49a7e14dcfSSatish Balay   PetscReal      *x, *f, *l, *u, *fb;
50a7e14dcfSSatish Balay   PetscReal      xval, fval, lval, uval;
51a7e14dcfSSatish Balay   PetscErrorCode ierr;
52a7e14dcfSSatish Balay   PetscInt       low[5], high[5], n, i;
53a7e14dcfSSatish Balay 
54a7e14dcfSSatish Balay   PetscFunctionBegin;
55a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,1);
56a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID,2);
57a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID,3);
58a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID,4);
59a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID,4);
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr);
62a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr);
63a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr);
64a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr);
65a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr);
66a7e14dcfSSatish Balay 
67a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
682d0e5244SBarry Smith     if (low[0] != low[i] || high[0] != high[i]) SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
69a7e14dcfSSatish Balay   }
70a7e14dcfSSatish Balay 
71a7e14dcfSSatish Balay   ierr = VecGetArray(X, &x);CHKERRQ(ierr);
72a7e14dcfSSatish Balay   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
73a7e14dcfSSatish Balay   ierr = VecGetArray(L, &l);CHKERRQ(ierr);
74a7e14dcfSSatish Balay   ierr = VecGetArray(U, &u);CHKERRQ(ierr);
75a7e14dcfSSatish Balay   ierr = VecGetArray(FB, &fb);CHKERRQ(ierr);
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
80a7e14dcfSSatish Balay     xval = x[i]; fval = f[i];
81a7e14dcfSSatish Balay     lval = l[i]; uval = u[i];
82a7e14dcfSSatish Balay 
83e270355aSBarry Smith     if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) {
84a7e14dcfSSatish Balay       fb[i] = -fval;
85e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
86a7e14dcfSSatish Balay       fb[i] = -Fischer(uval - xval, -fval);
87e270355aSBarry Smith     } else if (uval >=  PETSC_INFINITY) {
88a7e14dcfSSatish Balay       fb[i] =  Fischer(xval - lval,  fval);
892d0e5244SBarry Smith     } else if (lval == uval) {
90a7e14dcfSSatish Balay       fb[i] = lval - xval;
912d0e5244SBarry Smith     } else {
92a7e14dcfSSatish Balay       fval  =  Fischer(uval - xval, -fval);
93a7e14dcfSSatish Balay       fb[i] =  Fischer(xval - lval,  fval);
94a7e14dcfSSatish Balay     }
95a7e14dcfSSatish Balay   }
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay   ierr = VecRestoreArray(X, &x);CHKERRQ(ierr);
98a7e14dcfSSatish Balay   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
99a7e14dcfSSatish Balay   ierr = VecRestoreArray(L, &l);CHKERRQ(ierr);
100a7e14dcfSSatish Balay   ierr = VecRestoreArray(U, &u);CHKERRQ(ierr);
101a7e14dcfSSatish Balay   ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr);
102a7e14dcfSSatish Balay   PetscFunctionReturn(0);
103a7e14dcfSSatish Balay }
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c)
106a7e14dcfSSatish Balay {
107a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
108a7e14dcfSSatish Balay    if (a + b <= 0) {
109a7e14dcfSSatish Balay      return PetscSqrtScalar(a*a + b*b + 2.0*c*c) - (a + b);
110a7e14dcfSSatish Balay    }
111a7e14dcfSSatish Balay    return 2.0*(c*c - a*b) / (PetscSqrtScalar(a*a + b*b + 2.0*c*c) + (a + b));
112a7e14dcfSSatish Balay }
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay #undef __FUNCT__
115a7e14dcfSSatish Balay #define __FUNCT__ "VecSFischer"
116a7e14dcfSSatish Balay /*@
117a7e14dcfSSatish Balay    VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for
118a7e14dcfSSatish Balay    complementarity problems.
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay    Logically Collective on vectors
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay    Input Parameters:
123a7e14dcfSSatish Balay +  X - current point
124a7e14dcfSSatish Balay .  F - function evaluated at x
125a7e14dcfSSatish Balay .  L - lower bounds
126a7e14dcfSSatish Balay .  U - upper bounds
127a7e14dcfSSatish Balay -  mu - smoothing parameter
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay    Output Parameters:
130a7e14dcfSSatish Balay .  FB - The Smoothed Fischer-Burmeister function vector
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay    Notes:
133a7e14dcfSSatish Balay    The Smoothed Fischer-Burmeister function is defined as
134a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b
135a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
136a7e14dcfSSatish Balay    system of equations.
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay    The result of this function is done by cases:
139a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i] - 2*mu*x[i]
140a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i], mu)
141a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i], mu)
142a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu)
143a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay    Level: developer
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay .seealso  VecFischer()
148a7e14dcfSSatish Balay @*/
149a7e14dcfSSatish Balay PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB)
150a7e14dcfSSatish Balay {
151a7e14dcfSSatish Balay   PetscReal      *x, *f, *l, *u, *fb;
152a7e14dcfSSatish Balay   PetscReal      xval, fval, lval, uval;
153a7e14dcfSSatish Balay   PetscErrorCode ierr;
154a7e14dcfSSatish Balay   PetscInt       low[5], high[5], n, i;
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay   PetscFunctionBegin;
157a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,1);
158a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID,2);
159a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID,3);
160a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID,4);
161a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID,6);
162a7e14dcfSSatish Balay 
163a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr);
164a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr);
165a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr);
166a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr);
167a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr);
168a7e14dcfSSatish Balay 
169a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
1702d0e5244SBarry Smith     if (low[0] != low[i] || high[0] != high[i]) SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
171a7e14dcfSSatish Balay   }
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay   ierr = VecGetArray(X, &x);CHKERRQ(ierr);
174a7e14dcfSSatish Balay   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
175a7e14dcfSSatish Balay   ierr = VecGetArray(L, &l);CHKERRQ(ierr);
176a7e14dcfSSatish Balay   ierr = VecGetArray(U, &u);CHKERRQ(ierr);
177a7e14dcfSSatish Balay   ierr = VecGetArray(FB, &fb);CHKERRQ(ierr);
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
182a7e14dcfSSatish Balay     xval = (*x++); fval = (*f++);
183a7e14dcfSSatish Balay     lval = (*l++); uval = (*u++);
184a7e14dcfSSatish Balay 
185e270355aSBarry Smith     if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) {
186a7e14dcfSSatish Balay       (*fb++) = -fval - mu*xval;
187e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
188a7e14dcfSSatish Balay       (*fb++) = -SFischer(uval - xval, -fval, mu);
189e270355aSBarry Smith     } else if (uval >=  PETSC_INFINITY) {
190a7e14dcfSSatish Balay       (*fb++) =  SFischer(xval - lval,  fval, mu);
1912d0e5244SBarry Smith     } else if (lval == uval) {
192a7e14dcfSSatish Balay       (*fb++) = lval - xval;
1932d0e5244SBarry Smith     } else {
194a7e14dcfSSatish Balay       fval    =  SFischer(uval - xval, -fval, mu);
195a7e14dcfSSatish Balay       (*fb++) =  SFischer(xval - lval,  fval, mu);
196a7e14dcfSSatish Balay     }
197a7e14dcfSSatish Balay   }
198a7e14dcfSSatish Balay   x -= n; f -= n; l -=n; u -= n; fb -= n;
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay   ierr = VecRestoreArray(X, &x);CHKERRQ(ierr);
201a7e14dcfSSatish Balay   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
202a7e14dcfSSatish Balay   ierr = VecRestoreArray(L, &l);CHKERRQ(ierr);
203a7e14dcfSSatish Balay   ierr = VecRestoreArray(U, &u);CHKERRQ(ierr);
204a7e14dcfSSatish Balay   ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr);
205a7e14dcfSSatish Balay   PetscFunctionReturn(0);
206a7e14dcfSSatish Balay }
207a7e14dcfSSatish Balay 
208a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal fischnorm(PetscReal a, PetscReal b)
209a7e14dcfSSatish Balay {
210a7e14dcfSSatish Balay   return PetscSqrtScalar(a*a + b*b);
211a7e14dcfSSatish Balay }
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c)
214a7e14dcfSSatish Balay {
215a7e14dcfSSatish Balay   return PetscSqrtScalar(a*a + b*b + 2.0*c*c);
216a7e14dcfSSatish Balay }
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay #undef __FUNCT__
219a7e14dcfSSatish Balay #define __FUNCT__ "D_Fischer"
220a7e14dcfSSatish Balay /*@
221a7e14dcfSSatish Balay    D_Fischer - Calculates an element of the B-subdifferential of the
222a7e14dcfSSatish Balay    Fischer-Burmeister function for complementarity problems.
223a7e14dcfSSatish Balay 
224a7e14dcfSSatish Balay    Collective on jac
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay    Input Parameters:
227a7e14dcfSSatish Balay +  jac - the jacobian of f at X
228a7e14dcfSSatish Balay .  X - current point
229a7e14dcfSSatish Balay .  Con - constraints function evaluated at X
230a7e14dcfSSatish Balay .  XL - lower bounds
231a7e14dcfSSatish Balay .  XU - upper bounds
232a7e14dcfSSatish Balay .  t1 - work vector
233a7e14dcfSSatish Balay -  t2 - work vector
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay    Output Parameters:
236a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
237a7e14dcfSSatish Balay -  Db - row scaling component of the result
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay    Level: developer
240a7e14dcfSSatish Balay 
241a7e14dcfSSatish Balay .seealso: VecFischer()
242a7e14dcfSSatish Balay @*/
2432d0e5244SBarry Smith PetscErrorCode D_Fischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db)
244a7e14dcfSSatish Balay {
245a7e14dcfSSatish Balay   PetscErrorCode ierr;
246a7e14dcfSSatish Balay   PetscInt       i,nn;
247a7e14dcfSSatish Balay   PetscReal      *x,*f,*l,*u,*da,*db,*t1,*t2;
248a7e14dcfSSatish Balay   PetscReal      ai,bi,ci,di,ei;
249a7e14dcfSSatish Balay 
250a7e14dcfSSatish Balay   PetscFunctionBegin;
251a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr);
252a7e14dcfSSatish Balay   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
253a7e14dcfSSatish Balay   ierr = VecGetArray(Con,&f);CHKERRQ(ierr);
254a7e14dcfSSatish Balay   ierr = VecGetArray(XL,&l);CHKERRQ(ierr);
255a7e14dcfSSatish Balay   ierr = VecGetArray(XU,&u);CHKERRQ(ierr);
256a7e14dcfSSatish Balay   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
257a7e14dcfSSatish Balay   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
258a7e14dcfSSatish Balay   ierr = VecGetArray(T1,&t1);CHKERRQ(ierr);
259a7e14dcfSSatish Balay   ierr = VecGetArray(T2,&t2);CHKERRQ(ierr);
260a7e14dcfSSatish Balay 
261a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
262a7e14dcfSSatish Balay     da[i] = 0;
263a7e14dcfSSatish Balay     db[i] = 0;
264a7e14dcfSSatish Balay     t1[i] = 0;
265a7e14dcfSSatish Balay 
266a7e14dcfSSatish Balay     if (PetscAbsReal(f[i]) <= PETSC_MACHINE_EPSILON) {
267e270355aSBarry Smith       if (l[i] > PETSC_NINFINITY && PetscAbsReal(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) {
268a7e14dcfSSatish Balay         t1[i] = 1;
269a7e14dcfSSatish Balay         da[i] = 1;
270a7e14dcfSSatish Balay       }
271a7e14dcfSSatish Balay 
272e270355aSBarry Smith       if (u[i] <  PETSC_INFINITY && PetscAbsReal(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) {
273a7e14dcfSSatish Balay         t1[i] = 1;
274a7e14dcfSSatish Balay         db[i] = 1;
275a7e14dcfSSatish Balay       }
276a7e14dcfSSatish Balay     }
277a7e14dcfSSatish Balay   }
278a7e14dcfSSatish Balay 
279a7e14dcfSSatish Balay   ierr = VecRestoreArray(T1,&t1);CHKERRQ(ierr);
280a7e14dcfSSatish Balay   ierr = VecRestoreArray(T2,&t2);CHKERRQ(ierr);
281a7e14dcfSSatish Balay   ierr = MatMult(jac,T1,T2);CHKERRQ(ierr);
282a7e14dcfSSatish Balay   ierr = VecGetArray(T2,&t2);CHKERRQ(ierr);
283a7e14dcfSSatish Balay 
284a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
285e270355aSBarry Smith     if ((l[i] <= PETSC_NINFINITY) && (u[i] >= PETSC_INFINITY)) {
286a7e14dcfSSatish Balay       da[i] = 0;
287a7e14dcfSSatish Balay       db[i] = -1;
288e270355aSBarry Smith     } else if (l[i] <= PETSC_NINFINITY) {
289a7e14dcfSSatish Balay       if (db[i] >= 1) {
290a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
291a7e14dcfSSatish Balay 
292a7e14dcfSSatish Balay         da[i] = -1/ai - 1;
293a7e14dcfSSatish Balay         db[i] = -t2[i]/ai - 1;
2942d0e5244SBarry Smith       } else {
295a7e14dcfSSatish Balay         bi = u[i] - x[i];
296a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
297a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
298a7e14dcfSSatish Balay 
299a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
300a7e14dcfSSatish Balay         db[i] = -f[i] / ai - 1;
301a7e14dcfSSatish Balay       }
302e270355aSBarry Smith     } else if (u[i] >=  PETSC_INFINITY) {
303a7e14dcfSSatish Balay       if (da[i] >= 1) {
304a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay         da[i] = 1 / ai - 1;
307a7e14dcfSSatish Balay         db[i] = t2[i] / ai - 1;
3082d0e5244SBarry Smith       } else {
309a7e14dcfSSatish Balay         bi = x[i] - l[i];
310a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
311a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
312a7e14dcfSSatish Balay 
313a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
314a7e14dcfSSatish Balay         db[i] = f[i] / ai - 1;
315a7e14dcfSSatish Balay       }
3162d0e5244SBarry Smith     } else if (l[i] == u[i]) {
317a7e14dcfSSatish Balay       da[i] = -1;
318a7e14dcfSSatish Balay       db[i] = 0;
3192d0e5244SBarry Smith     } else {
320a7e14dcfSSatish Balay       if (db[i] >= 1) {
321a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
322a7e14dcfSSatish Balay 
323a7e14dcfSSatish Balay         ci = 1 / ai + 1;
324a7e14dcfSSatish Balay         di = t2[i] / ai + 1;
3252d0e5244SBarry Smith       } else {
326a7e14dcfSSatish Balay         bi = x[i] - u[i];
327a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
328a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
329a7e14dcfSSatish Balay 
330a7e14dcfSSatish Balay         ci = bi / ai + 1;
331a7e14dcfSSatish Balay         di = f[i] / ai + 1;
332a7e14dcfSSatish Balay       }
333a7e14dcfSSatish Balay 
334a7e14dcfSSatish Balay       if (da[i] >= 1) {
335a7e14dcfSSatish Balay         bi = ci + di*t2[i];
336a7e14dcfSSatish Balay         ai = fischnorm(1, bi);
337a7e14dcfSSatish Balay 
338a7e14dcfSSatish Balay         bi = bi / ai - 1;
339a7e14dcfSSatish Balay         ai = 1 / ai - 1;
3402d0e5244SBarry Smith       } else {
341a7e14dcfSSatish Balay         ei = Fischer(u[i] - x[i], -f[i]);
342a7e14dcfSSatish Balay         ai = fischnorm(x[i] - l[i], ei);
343a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
344a7e14dcfSSatish Balay 
345a7e14dcfSSatish Balay         bi = ei / ai - 1;
346a7e14dcfSSatish Balay         ai = (x[i] - l[i]) / ai - 1;
347a7e14dcfSSatish Balay       }
348a7e14dcfSSatish Balay 
349a7e14dcfSSatish Balay       da[i] = ai + bi*ci;
350a7e14dcfSSatish Balay       db[i] = bi*di;
351a7e14dcfSSatish Balay     }
352a7e14dcfSSatish Balay   }
353a7e14dcfSSatish Balay 
354a7e14dcfSSatish Balay   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
355a7e14dcfSSatish Balay   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
356a7e14dcfSSatish Balay   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
357a7e14dcfSSatish Balay   ierr = VecRestoreArray(Con,&f);CHKERRQ(ierr);
358a7e14dcfSSatish Balay   ierr = VecRestoreArray(XL,&l);CHKERRQ(ierr);
359a7e14dcfSSatish Balay   ierr = VecRestoreArray(XU,&u);CHKERRQ(ierr);
360a7e14dcfSSatish Balay   ierr = VecRestoreArray(T2,&t2);CHKERRQ(ierr);
361a7e14dcfSSatish Balay   PetscFunctionReturn(0);
3628e3154b5SSatish Balay }
363a7e14dcfSSatish Balay 
364a7e14dcfSSatish Balay #undef __FUNCT__
365a7e14dcfSSatish Balay #define __FUNCT__ "D_SFischer"
366a7e14dcfSSatish Balay /*@
367a7e14dcfSSatish Balay    D_SFischer - Calculates an element of the B-subdifferential of the
368a7e14dcfSSatish Balay    smoothed Fischer-Burmeister function for complementarity problems.
369a7e14dcfSSatish Balay 
370a7e14dcfSSatish Balay    Collective on jac
371a7e14dcfSSatish Balay 
372a7e14dcfSSatish Balay    Input Parameters:
373a7e14dcfSSatish Balay +  jac - the jacobian of f at X
374a7e14dcfSSatish Balay .  X - current point
375a7e14dcfSSatish Balay .  F - constraint function evaluated at X
376a7e14dcfSSatish Balay .  XL - lower bounds
377a7e14dcfSSatish Balay .  XU - upper bounds
378a7e14dcfSSatish Balay .  mu - smoothing parameter
379a7e14dcfSSatish Balay .  T1 - work vector
380a7e14dcfSSatish Balay -  T2 - work vector
381a7e14dcfSSatish Balay 
382a7e14dcfSSatish Balay    Output Parameter:
383a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
384a7e14dcfSSatish Balay .  Db - row scaling component of the result
385a7e14dcfSSatish Balay -  Dm - derivative with respect to scaling parameter
386a7e14dcfSSatish Balay 
387a7e14dcfSSatish Balay    Level: developer
388a7e14dcfSSatish Balay 
389a7e14dcfSSatish Balay .seealso D_Fischer()
390a7e14dcfSSatish Balay @*/
3912d0e5244SBarry Smith PetscErrorCode D_SFischer(Mat jac, Vec X, Vec Con,Vec XL, Vec XU, PetscReal mu,Vec T1, Vec T2,Vec Da, Vec Db, Vec Dm)
392a7e14dcfSSatish Balay {
393a7e14dcfSSatish Balay   PetscErrorCode ierr;
394a7e14dcfSSatish Balay   PetscInt       i,nn;
395a7e14dcfSSatish Balay   PetscReal      *x, *f, *l, *u, *da, *db, *dm;
396a7e14dcfSSatish Balay   PetscReal      ai, bi, ci, di, ei, fi;
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay   PetscFunctionBegin;
399a7e14dcfSSatish Balay   if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) {
400a7e14dcfSSatish Balay     ierr = VecZeroEntries(Dm);CHKERRQ(ierr);
401a7e14dcfSSatish Balay     ierr = D_Fischer(jac, X, Con, XL, XU, T1, T2, Da, Db);CHKERRQ(ierr);
4022d0e5244SBarry Smith   } else {
403a7e14dcfSSatish Balay     ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr);
404a7e14dcfSSatish Balay     ierr = VecGetArray(X,&x);CHKERRQ(ierr);
405a7e14dcfSSatish Balay     ierr = VecGetArray(Con,&f);CHKERRQ(ierr);
406a7e14dcfSSatish Balay     ierr = VecGetArray(XL,&l);CHKERRQ(ierr);
407a7e14dcfSSatish Balay     ierr = VecGetArray(XU,&u);CHKERRQ(ierr);
408a7e14dcfSSatish Balay     ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
409a7e14dcfSSatish Balay     ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
410a7e14dcfSSatish Balay     ierr = VecGetArray(Dm,&dm);CHKERRQ(ierr);
411a7e14dcfSSatish Balay 
412a7e14dcfSSatish Balay     for (i = 0; i < nn; ++i) {
413e270355aSBarry Smith       if ((l[i] <= PETSC_NINFINITY) && (u[i] >= PETSC_INFINITY)) {
414a7e14dcfSSatish Balay         da[i] = -mu;
415a7e14dcfSSatish Balay         db[i] = -1;
416a7e14dcfSSatish Balay         dm[i] = -x[i];
417e270355aSBarry Smith       } else if (l[i] <= PETSC_NINFINITY) {
418a7e14dcfSSatish Balay         bi = u[i] - x[i];
419a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
420a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
421a7e14dcfSSatish Balay 
422a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
423a7e14dcfSSatish Balay         db[i] = -f[i] / ai - 1;
424a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
425e270355aSBarry Smith       } else if (u[i] >=  PETSC_INFINITY) {
426a7e14dcfSSatish Balay         bi = x[i] - l[i];
427a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
428a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
431a7e14dcfSSatish Balay         db[i] = f[i] / ai - 1;
432a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
4332d0e5244SBarry Smith       } else if (l[i] == u[i]) {
434a7e14dcfSSatish Balay         da[i] = -1;
435a7e14dcfSSatish Balay         db[i] = 0;
436a7e14dcfSSatish Balay         dm[i] = 0;
4372d0e5244SBarry Smith       } else {
438a7e14dcfSSatish Balay         bi = x[i] - u[i];
439a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
440a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
441a7e14dcfSSatish Balay 
442a7e14dcfSSatish Balay         ci = bi / ai + 1;
443a7e14dcfSSatish Balay         di = f[i] / ai + 1;
444a7e14dcfSSatish Balay         fi = 2.0 * mu / ai;
445a7e14dcfSSatish Balay 
446a7e14dcfSSatish Balay         ei = SFischer(u[i] - x[i], -f[i], mu);
447a7e14dcfSSatish Balay         ai = fischsnorm(x[i] - l[i], ei, mu);
448a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
449a7e14dcfSSatish Balay 
450a7e14dcfSSatish Balay         bi = ei / ai - 1;
451a7e14dcfSSatish Balay         ei = 2.0 * mu / ei;
452a7e14dcfSSatish Balay         ai = (x[i] - l[i]) / ai - 1;
453a7e14dcfSSatish Balay 
454a7e14dcfSSatish Balay         da[i] = ai + bi*ci;
455a7e14dcfSSatish Balay         db[i] = bi*di;
456a7e14dcfSSatish Balay         dm[i] = ei + bi*fi;
457a7e14dcfSSatish Balay       }
458a7e14dcfSSatish Balay     }
459a7e14dcfSSatish Balay 
460a7e14dcfSSatish Balay     ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
461a7e14dcfSSatish Balay     ierr = VecRestoreArray(Con,&f);CHKERRQ(ierr);
462a7e14dcfSSatish Balay     ierr = VecRestoreArray(XL,&l);CHKERRQ(ierr);
463a7e14dcfSSatish Balay     ierr = VecRestoreArray(XU,&u);CHKERRQ(ierr);
464a7e14dcfSSatish Balay     ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
465a7e14dcfSSatish Balay     ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
466a7e14dcfSSatish Balay     ierr = VecRestoreArray(Dm,&dm);CHKERRQ(ierr);
467a7e14dcfSSatish Balay   }
468a7e14dcfSSatish Balay   PetscFunctionReturn(0);
469a7e14dcfSSatish Balay }
470a7e14dcfSSatish Balay 
471