xref: /petsc/src/tao/util/tao_util.c (revision 46bdf8c8e7980dd1dca8bd7d868b63fbdf2185ad)
1aaa7dc30SBarry Smith #include <petsc-private/petscimpl.h>
2ba92ff59SBarry 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) {
9*46bdf8c8SLisandro Dalcin      return PetscSqrtReal(a*a + b*b) - (a + b);
10a7e14dcfSSatish Balay    }
11*46bdf8c8SLisandro Dalcin    return -2.0*a*b / (PetscSqrtReal(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 {
49*46bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
50*46bdf8c8SLisandro Dalcin   PetscScalar       *fb;
51a7e14dcfSSatish Balay   PetscReal         xval, fval, lval, uval;
52a7e14dcfSSatish Balay   PetscErrorCode    ierr;
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);
58a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID,3);
59a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID,4);
60a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID,4);
61a7e14dcfSSatish Balay 
62a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr);
63a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr);
64a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr);
65a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr);
66a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr);
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
692d0e5244SBarry Smith     if (low[0] != low[i] || high[0] != high[i]) SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
70a7e14dcfSSatish Balay   }
71a7e14dcfSSatish Balay 
725e081366SBarry Smith   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
735e081366SBarry Smith   ierr = VecGetArrayRead(F, &f);CHKERRQ(ierr);
745e081366SBarry Smith   ierr = VecGetArrayRead(L, &l);CHKERRQ(ierr);
755e081366SBarry Smith   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
76a7e14dcfSSatish Balay   ierr = VecGetArray(FB, &fb);CHKERRQ(ierr);
77a7e14dcfSSatish Balay 
78a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
81a7e14dcfSSatish Balay     xval = x[i]; fval = f[i];
82a7e14dcfSSatish Balay     lval = l[i]; uval = u[i];
83a7e14dcfSSatish Balay 
84e270355aSBarry Smith     if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) {
85a7e14dcfSSatish Balay       fb[i] = -fval;
86e270355aSBarry Smith     } else if (lval <= -PETSC_INFINITY) {
87a7e14dcfSSatish Balay       fb[i] = -Fischer(uval - xval, -fval);
88e270355aSBarry Smith     } else if (uval >=  PETSC_INFINITY) {
89a7e14dcfSSatish Balay       fb[i] =  Fischer(xval - lval,  fval);
902d0e5244SBarry Smith     } else if (lval == uval) {
91a7e14dcfSSatish Balay       fb[i] = lval - xval;
922d0e5244SBarry Smith     } else {
93a7e14dcfSSatish Balay       fval  =  Fischer(uval - xval, -fval);
94a7e14dcfSSatish Balay       fb[i] =  Fischer(xval - lval,  fval);
95a7e14dcfSSatish Balay     }
96a7e14dcfSSatish Balay   }
97a7e14dcfSSatish Balay 
985e081366SBarry Smith   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
995e081366SBarry Smith   ierr = VecRestoreArrayRead(F, &f);CHKERRQ(ierr);
1005e081366SBarry Smith   ierr = VecRestoreArrayRead(L, &l);CHKERRQ(ierr);
1015e081366SBarry Smith   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
102a7e14dcfSSatish Balay   ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr);
103a7e14dcfSSatish Balay   PetscFunctionReturn(0);
104a7e14dcfSSatish Balay }
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c)
107a7e14dcfSSatish Balay {
108a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
109a7e14dcfSSatish Balay    if (a + b <= 0) {
110a7e14dcfSSatish Balay      return PetscSqrtScalar(a*a + b*b + 2.0*c*c) - (a + b);
111a7e14dcfSSatish Balay    }
112a7e14dcfSSatish Balay    return 2.0*(c*c - a*b) / (PetscSqrtScalar(a*a + b*b + 2.0*c*c) + (a + b));
113a7e14dcfSSatish Balay }
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay #undef __FUNCT__
116a7e14dcfSSatish Balay #define __FUNCT__ "VecSFischer"
117a7e14dcfSSatish Balay /*@
118a7e14dcfSSatish Balay    VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for
119a7e14dcfSSatish Balay    complementarity problems.
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay    Logically Collective on vectors
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay    Input Parameters:
124a7e14dcfSSatish Balay +  X - current point
125a7e14dcfSSatish Balay .  F - function evaluated at x
126a7e14dcfSSatish Balay .  L - lower bounds
127a7e14dcfSSatish Balay .  U - upper bounds
128a7e14dcfSSatish Balay -  mu - smoothing parameter
129a7e14dcfSSatish Balay 
130a7e14dcfSSatish Balay    Output Parameters:
131a7e14dcfSSatish Balay .  FB - The Smoothed Fischer-Burmeister function vector
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay    Notes:
134a7e14dcfSSatish Balay    The Smoothed Fischer-Burmeister function is defined as
135a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b
136a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
137a7e14dcfSSatish Balay    system of equations.
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay    The result of this function is done by cases:
140a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i] - 2*mu*x[i]
141a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i], mu)
142a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i], mu)
143a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu)
144a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
145a7e14dcfSSatish Balay 
146a7e14dcfSSatish Balay    Level: developer
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay .seealso  VecFischer()
149a7e14dcfSSatish Balay @*/
150a7e14dcfSSatish Balay PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB)
151a7e14dcfSSatish Balay {
152*46bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
153*46bdf8c8SLisandro Dalcin   PetscScalar       *fb;
154a7e14dcfSSatish Balay   PetscReal         xval, fval, lval, uval;
155a7e14dcfSSatish Balay   PetscErrorCode    ierr;
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 
165a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr);
166a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr);
167a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr);
168a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr);
169a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr);
170a7e14dcfSSatish Balay 
171a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
1722d0e5244SBarry Smith     if (low[0] != low[i] || high[0] != high[i]) SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
173a7e14dcfSSatish Balay   }
174a7e14dcfSSatish Balay 
1755e081366SBarry Smith   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
1765e081366SBarry Smith   ierr = VecGetArrayRead(F, &f);CHKERRQ(ierr);
1775e081366SBarry Smith   ierr = VecGetArrayRead(L, &l);CHKERRQ(ierr);
1785e081366SBarry Smith   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
179a7e14dcfSSatish Balay   ierr = VecGetArray(FB, &fb);CHKERRQ(ierr);
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
184a7e14dcfSSatish Balay     xval = (*x++); fval = (*f++);
185a7e14dcfSSatish Balay     lval = (*l++); uval = (*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   }
200a7e14dcfSSatish Balay   x -= n; f -= n; l -=n; u -= n; fb -= n;
201a7e14dcfSSatish Balay 
2025e081366SBarry Smith   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
2035e081366SBarry Smith   ierr = VecRestoreArrayRead(F, &f);CHKERRQ(ierr);
2045e081366SBarry Smith   ierr = VecRestoreArrayRead(L, &l);CHKERRQ(ierr);
2055e081366SBarry Smith   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
206a7e14dcfSSatish Balay   ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr);
207a7e14dcfSSatish Balay   PetscFunctionReturn(0);
208a7e14dcfSSatish Balay }
209a7e14dcfSSatish Balay 
210a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal fischnorm(PetscReal a, PetscReal b)
211a7e14dcfSSatish Balay {
212a7e14dcfSSatish Balay   return PetscSqrtScalar(a*a + b*b);
213a7e14dcfSSatish Balay }
214a7e14dcfSSatish Balay 
215a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c)
216a7e14dcfSSatish Balay {
217a7e14dcfSSatish Balay   return PetscSqrtScalar(a*a + b*b + 2.0*c*c);
218a7e14dcfSSatish Balay }
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay #undef __FUNCT__
221235fd6e6SBarry Smith #define __FUNCT__ "MatDFischer"
222a7e14dcfSSatish Balay /*@
223235fd6e6SBarry Smith    MatDFischer - Calculates an element of the B-subdifferential of the
224a7e14dcfSSatish Balay    Fischer-Burmeister function for complementarity problems.
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay    Collective on jac
227a7e14dcfSSatish Balay 
228a7e14dcfSSatish Balay    Input Parameters:
229a7e14dcfSSatish Balay +  jac - the jacobian of f at X
230a7e14dcfSSatish Balay .  X - current point
231a7e14dcfSSatish Balay .  Con - constraints function evaluated at X
232a7e14dcfSSatish Balay .  XL - lower bounds
233a7e14dcfSSatish Balay .  XU - upper bounds
234a7e14dcfSSatish Balay .  t1 - work vector
235a7e14dcfSSatish Balay -  t2 - work vector
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay    Output Parameters:
238a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
239a7e14dcfSSatish Balay -  Db - row scaling component of the result
240a7e14dcfSSatish Balay 
241a7e14dcfSSatish Balay    Level: developer
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay .seealso: VecFischer()
244a7e14dcfSSatish Balay @*/
245235fd6e6SBarry Smith PetscErrorCode MatDFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db)
246a7e14dcfSSatish Balay {
247a7e14dcfSSatish Balay   PetscErrorCode    ierr;
248a7e14dcfSSatish Balay   PetscInt          i,nn;
249*46bdf8c8SLisandro Dalcin   const PetscScalar *x,*f,*l,*u,*t2;
250*46bdf8c8SLisandro Dalcin   PetscScalar       *da,*db,*t1;
251a7e14dcfSSatish Balay   PetscReal          ai,bi,ci,di,ei;
252a7e14dcfSSatish Balay 
253a7e14dcfSSatish Balay   PetscFunctionBegin;
254a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr);
2555e081366SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
2565e081366SBarry Smith   ierr = VecGetArrayRead(Con,&f);CHKERRQ(ierr);
2575e081366SBarry Smith   ierr = VecGetArrayRead(XL,&l);CHKERRQ(ierr);
2585e081366SBarry Smith   ierr = VecGetArrayRead(XU,&u);CHKERRQ(ierr);
259a7e14dcfSSatish Balay   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
260a7e14dcfSSatish Balay   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
261a7e14dcfSSatish Balay   ierr = VecGetArray(T1,&t1);CHKERRQ(ierr);
2625e081366SBarry Smith   ierr = VecGetArrayRead(T2,&t2);CHKERRQ(ierr);
263a7e14dcfSSatish Balay 
264a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
265*46bdf8c8SLisandro Dalcin     da[i] = 0.0;
266*46bdf8c8SLisandro Dalcin     db[i] = 0.0;
267*46bdf8c8SLisandro Dalcin     t1[i] = 0.0;
268a7e14dcfSSatish Balay 
269*46bdf8c8SLisandro Dalcin     if (PetscAbsScalar(f[i]) <= PETSC_MACHINE_EPSILON) {
270*46bdf8c8SLisandro Dalcin       if (PetscRealPart(l[i]) > PETSC_NINFINITY && PetscAbsScalar(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) {
271*46bdf8c8SLisandro Dalcin         t1[i] = 1.0;
272*46bdf8c8SLisandro Dalcin         da[i] = 1.0;
273a7e14dcfSSatish Balay       }
274a7e14dcfSSatish Balay 
275*46bdf8c8SLisandro Dalcin       if (PetscRealPart(u[i]) <  PETSC_INFINITY && PetscAbsScalar(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) {
276*46bdf8c8SLisandro Dalcin         t1[i] = 1.0;
277*46bdf8c8SLisandro Dalcin         db[i] = 1.0;
278a7e14dcfSSatish Balay       }
279a7e14dcfSSatish Balay     }
280a7e14dcfSSatish Balay   }
281a7e14dcfSSatish Balay 
282a7e14dcfSSatish Balay   ierr = VecRestoreArray(T1,&t1);CHKERRQ(ierr);
2835e081366SBarry Smith   ierr = VecRestoreArrayRead(T2,&t2);CHKERRQ(ierr);
284a7e14dcfSSatish Balay   ierr = MatMult(jac,T1,T2);CHKERRQ(ierr);
2855e081366SBarry Smith   ierr = VecGetArrayRead(T2,&t2);CHKERRQ(ierr);
286a7e14dcfSSatish Balay 
287a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
288*46bdf8c8SLisandro Dalcin     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) {
289*46bdf8c8SLisandro Dalcin       da[i] = 0.0;
290*46bdf8c8SLisandro Dalcin       db[i] = -1.0;
291*46bdf8c8SLisandro Dalcin     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {
292*46bdf8c8SLisandro Dalcin       if (PetscRealPart(db[i]) >= 1) {
293a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
294a7e14dcfSSatish Balay 
295*46bdf8c8SLisandro Dalcin         da[i] = -1.0 / ai - 1.0;
296*46bdf8c8SLisandro Dalcin         db[i] = -t2[i] / ai - 1.0;
2972d0e5244SBarry Smith       } else {
298a7e14dcfSSatish Balay         bi = u[i] - x[i];
299a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
300a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
301a7e14dcfSSatish Balay 
302*46bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
303*46bdf8c8SLisandro Dalcin         db[i] = -f[i] / ai - 1.0;
304a7e14dcfSSatish Balay       }
305*46bdf8c8SLisandro Dalcin     } else if (PetscRealPart(u[i]) >=  PETSC_INFINITY) {
306*46bdf8c8SLisandro Dalcin       if (PetscRealPart(da[i]) >= 1) {
307a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
308a7e14dcfSSatish Balay 
309*46bdf8c8SLisandro Dalcin         da[i] = 1.0 / ai - 1.0;
310*46bdf8c8SLisandro Dalcin         db[i] = t2[i] / ai - 1.0;
3112d0e5244SBarry Smith       } else {
312a7e14dcfSSatish Balay         bi = x[i] - l[i];
313a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
314a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
315a7e14dcfSSatish Balay 
316*46bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
317*46bdf8c8SLisandro Dalcin         db[i] = f[i] / ai - 1.0;
318a7e14dcfSSatish Balay       }
3192d0e5244SBarry Smith     } else if (l[i] == u[i]) {
320*46bdf8c8SLisandro Dalcin       da[i] = -1.0;
321*46bdf8c8SLisandro Dalcin       db[i] = 0.0;
3222d0e5244SBarry Smith     } else {
323*46bdf8c8SLisandro Dalcin       if (PetscRealPart(db[i]) >= 1) {
324a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
325a7e14dcfSSatish Balay 
326*46bdf8c8SLisandro Dalcin         ci = 1.0 / ai + 1.0;
327*46bdf8c8SLisandro Dalcin         di = t2[i] / ai + 1.0;
3282d0e5244SBarry Smith       } else {
329a7e14dcfSSatish Balay         bi = x[i] - u[i];
330a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
331a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
332a7e14dcfSSatish Balay 
333*46bdf8c8SLisandro Dalcin         ci = bi / ai + 1.0;
334*46bdf8c8SLisandro Dalcin         di = f[i] / ai + 1.0;
335a7e14dcfSSatish Balay       }
336a7e14dcfSSatish Balay 
337*46bdf8c8SLisandro Dalcin       if (PetscRealPart(da[i]) >= 1) {
338a7e14dcfSSatish Balay         bi = ci + di*t2[i];
339a7e14dcfSSatish Balay         ai = fischnorm(1, bi);
340a7e14dcfSSatish Balay 
341*46bdf8c8SLisandro Dalcin         bi = bi / ai - 1.0;
342*46bdf8c8SLisandro Dalcin         ai = 1.0 / ai - 1.0;
3432d0e5244SBarry Smith       } else {
344a7e14dcfSSatish Balay         ei = Fischer(u[i] - x[i], -f[i]);
345a7e14dcfSSatish Balay         ai = fischnorm(x[i] - l[i], ei);
346a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
347a7e14dcfSSatish Balay 
348*46bdf8c8SLisandro Dalcin         bi = ei / ai - 1.0;
349*46bdf8c8SLisandro Dalcin         ai = (x[i] - l[i]) / ai - 1.0;
350a7e14dcfSSatish Balay       }
351a7e14dcfSSatish Balay 
352a7e14dcfSSatish Balay       da[i] = ai + bi*ci;
353a7e14dcfSSatish Balay       db[i] = bi*di;
354a7e14dcfSSatish Balay     }
355a7e14dcfSSatish Balay   }
356a7e14dcfSSatish Balay 
357a7e14dcfSSatish Balay   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
358a7e14dcfSSatish Balay   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
3595e081366SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
3605e081366SBarry Smith   ierr = VecRestoreArrayRead(Con,&f);CHKERRQ(ierr);
3615e081366SBarry Smith   ierr = VecRestoreArrayRead(XL,&l);CHKERRQ(ierr);
3625e081366SBarry Smith   ierr = VecRestoreArrayRead(XU,&u);CHKERRQ(ierr);
3635e081366SBarry Smith   ierr = VecRestoreArrayRead(T2,&t2);CHKERRQ(ierr);
364a7e14dcfSSatish Balay   PetscFunctionReturn(0);
3658e3154b5SSatish Balay }
366a7e14dcfSSatish Balay 
367a7e14dcfSSatish Balay #undef __FUNCT__
368235fd6e6SBarry Smith #define __FUNCT__ "MatDSFischer"
369a7e14dcfSSatish Balay /*@
370235fd6e6SBarry Smith    MatDSFischer - Calculates an element of the B-subdifferential of the
371a7e14dcfSSatish Balay    smoothed Fischer-Burmeister function for complementarity problems.
372a7e14dcfSSatish Balay 
373a7e14dcfSSatish Balay    Collective on jac
374a7e14dcfSSatish Balay 
375a7e14dcfSSatish Balay    Input Parameters:
376a7e14dcfSSatish Balay +  jac - the jacobian of f at X
377a7e14dcfSSatish Balay .  X - current point
378a7e14dcfSSatish Balay .  F - constraint function evaluated at X
379a7e14dcfSSatish Balay .  XL - lower bounds
380a7e14dcfSSatish Balay .  XU - upper bounds
381a7e14dcfSSatish Balay .  mu - smoothing parameter
382a7e14dcfSSatish Balay .  T1 - work vector
383a7e14dcfSSatish Balay -  T2 - work vector
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay    Output Parameter:
386a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
387a7e14dcfSSatish Balay .  Db - row scaling component of the result
388a7e14dcfSSatish Balay -  Dm - derivative with respect to scaling parameter
389a7e14dcfSSatish Balay 
390a7e14dcfSSatish Balay    Level: developer
391a7e14dcfSSatish Balay 
392235fd6e6SBarry Smith .seealso MatDFischer()
393a7e14dcfSSatish Balay @*/
394235fd6e6SBarry 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)
395a7e14dcfSSatish Balay {
396a7e14dcfSSatish Balay   PetscErrorCode    ierr;
397a7e14dcfSSatish Balay   PetscInt          i,nn;
398*46bdf8c8SLisandro Dalcin   const PetscScalar *x, *f, *l, *u;
399*46bdf8c8SLisandro Dalcin   PetscScalar       *da, *db, *dm;
400a7e14dcfSSatish Balay   PetscReal         ai, bi, ci, di, ei, fi;
401a7e14dcfSSatish Balay 
402a7e14dcfSSatish Balay   PetscFunctionBegin;
403a7e14dcfSSatish Balay   if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) {
404a7e14dcfSSatish Balay     ierr = VecZeroEntries(Dm);CHKERRQ(ierr);
405235fd6e6SBarry Smith     ierr = MatDFischer(jac, X, Con, XL, XU, T1, T2, Da, Db);CHKERRQ(ierr);
4062d0e5244SBarry Smith   } else {
407a7e14dcfSSatish Balay     ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr);
4085e081366SBarry Smith     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
4095e081366SBarry Smith     ierr = VecGetArrayRead(Con,&f);CHKERRQ(ierr);
4105e081366SBarry Smith     ierr = VecGetArrayRead(XL,&l);CHKERRQ(ierr);
4115e081366SBarry Smith     ierr = VecGetArrayRead(XU,&u);CHKERRQ(ierr);
412a7e14dcfSSatish Balay     ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
413a7e14dcfSSatish Balay     ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
414a7e14dcfSSatish Balay     ierr = VecGetArray(Dm,&dm);CHKERRQ(ierr);
415a7e14dcfSSatish Balay 
416a7e14dcfSSatish Balay     for (i = 0; i < nn; ++i) {
417*46bdf8c8SLisandro Dalcin       if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) {
418a7e14dcfSSatish Balay         da[i] = -mu;
419*46bdf8c8SLisandro Dalcin         db[i] = -1.0;
420a7e14dcfSSatish Balay         dm[i] = -x[i];
421*46bdf8c8SLisandro Dalcin       } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {
422a7e14dcfSSatish Balay         bi = u[i] - x[i];
423a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
424a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
425a7e14dcfSSatish Balay 
426*46bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
427*46bdf8c8SLisandro Dalcin         db[i] = -f[i] / ai - 1.0;
428a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
429*46bdf8c8SLisandro Dalcin       } else if (PetscRealPart(u[i]) >=  PETSC_INFINITY) {
430a7e14dcfSSatish Balay         bi = x[i] - l[i];
431a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
432a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
433a7e14dcfSSatish Balay 
434*46bdf8c8SLisandro Dalcin         da[i] = bi / ai - 1.0;
435*46bdf8c8SLisandro Dalcin         db[i] = f[i] / ai - 1.0;
436a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
4372d0e5244SBarry Smith       } else if (l[i] == u[i]) {
438*46bdf8c8SLisandro Dalcin         da[i] = -1.0;
439*46bdf8c8SLisandro Dalcin         db[i] = 0.0;
440*46bdf8c8SLisandro Dalcin         dm[i] = 0.0;
4412d0e5244SBarry Smith       } else {
442a7e14dcfSSatish Balay         bi = x[i] - u[i];
443a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
444a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
445a7e14dcfSSatish Balay 
446*46bdf8c8SLisandro Dalcin         ci = bi / ai + 1.0;
447*46bdf8c8SLisandro Dalcin         di = f[i] / ai + 1.0;
448a7e14dcfSSatish Balay         fi = 2.0 * mu / ai;
449a7e14dcfSSatish Balay 
450a7e14dcfSSatish Balay         ei = SFischer(u[i] - x[i], -f[i], mu);
451a7e14dcfSSatish Balay         ai = fischsnorm(x[i] - l[i], ei, mu);
452a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
453a7e14dcfSSatish Balay 
454*46bdf8c8SLisandro Dalcin         bi = ei / ai - 1.0;
455a7e14dcfSSatish Balay         ei = 2.0 * mu / ei;
456*46bdf8c8SLisandro Dalcin         ai = (x[i] - l[i]) / ai - 1.0;
457a7e14dcfSSatish Balay 
458a7e14dcfSSatish Balay         da[i] = ai + bi*ci;
459a7e14dcfSSatish Balay         db[i] = bi*di;
460a7e14dcfSSatish Balay         dm[i] = ei + bi*fi;
461a7e14dcfSSatish Balay       }
462a7e14dcfSSatish Balay     }
463a7e14dcfSSatish Balay 
4645e081366SBarry Smith     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
4655e081366SBarry Smith     ierr = VecRestoreArrayRead(Con,&f);CHKERRQ(ierr);
4665e081366SBarry Smith     ierr = VecRestoreArrayRead(XL,&l);CHKERRQ(ierr);
4675e081366SBarry Smith     ierr = VecRestoreArrayRead(XU,&u);CHKERRQ(ierr);
468a7e14dcfSSatish Balay     ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
469a7e14dcfSSatish Balay     ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
470a7e14dcfSSatish Balay     ierr = VecRestoreArray(Dm,&dm);CHKERRQ(ierr);
471a7e14dcfSSatish Balay   }
472a7e14dcfSSatish Balay   PetscFunctionReturn(0);
473a7e14dcfSSatish Balay }
474a7e14dcfSSatish Balay 
475