xref: /petsc/src/tao/util/tao_util.c (revision aaa7dc30da3270cff6cb10b1db605b2ca746f216)
1*aaa7dc30SBarry Smith #include <petsc-private/petscimpl.h>
2*aaa7dc30SBarry Smith #include <tao.h> /*I "tao.h" I*/
3a7e14dcfSSatish Balay #include "tao_util.h" /*I "tao_util.h" I*/
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay #undef __FUNCT__
6a7e14dcfSSatish Balay #define __FUNCT__ "VecPow"
7a7e14dcfSSatish Balay /*@
8a7e14dcfSSatish Balay   VecPow - Replaces each component of a vector by x_i^p
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay   Logically Collective on v
11a7e14dcfSSatish Balay 
12a7e14dcfSSatish Balay   Input Parameter:
13a7e14dcfSSatish Balay + v - the vector
14a7e14dcfSSatish Balay - p - the exponent to use on each element
15a7e14dcfSSatish Balay 
16a7e14dcfSSatish Balay   Output Parameter:
17a7e14dcfSSatish Balay . v - the vector
18a7e14dcfSSatish Balay 
19a7e14dcfSSatish Balay   Level: intermediate
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay @*/
22a7e14dcfSSatish Balay PetscErrorCode VecPow(Vec v, PetscReal p)
23a7e14dcfSSatish Balay {
24a7e14dcfSSatish Balay   PetscErrorCode ierr;
25a7e14dcfSSatish Balay   PetscInt       n,i;
26a7e14dcfSSatish Balay   PetscReal      *v1;
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay   PetscFunctionBegin;
29a7e14dcfSSatish Balay   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay   ierr = VecGetArray(v, &v1);CHKERRQ(ierr);
32a7e14dcfSSatish Balay   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   if (1.0 == p) {
352d0e5244SBarry Smith   } else if (-1.0 == p) {
36a7e14dcfSSatish Balay     for (i = 0; i < n; ++i){
37a7e14dcfSSatish Balay       v1[i] = 1.0 / v1[i];
38a7e14dcfSSatish Balay     }
392d0e5244SBarry Smith   } else if (0.0 == p) {
40a7e14dcfSSatish Balay     for (i = 0; i < n; ++i){
41a7e14dcfSSatish Balay       /*  Not-a-number left alone
42a7e14dcfSSatish Balay           Infinity set to one  */
43a7e14dcfSSatish Balay       if (v1[i] == v1[i]) {
44a7e14dcfSSatish Balay         v1[i] = 1.0;
45a7e14dcfSSatish Balay       }
46a7e14dcfSSatish Balay     }
472d0e5244SBarry Smith   } else if (0.5 == p) {
48a7e14dcfSSatish Balay     for (i = 0; i < n; ++i) {
49a7e14dcfSSatish Balay       if (v1[i] >= 0) {
50a7e14dcfSSatish Balay         v1[i] = PetscSqrtScalar(v1[i]);
512d0e5244SBarry Smith       } else {
52a7e14dcfSSatish Balay         v1[i] = TAO_INFINITY;
53a7e14dcfSSatish Balay       }
54a7e14dcfSSatish Balay     }
552d0e5244SBarry Smith   } else if (-0.5 == p) {
56a7e14dcfSSatish Balay     for (i = 0; i < n; ++i) {
57a7e14dcfSSatish Balay       if (v1[i] >= 0) {
58a7e14dcfSSatish Balay         v1[i] = 1.0 / PetscSqrtScalar(v1[i]);
592d0e5244SBarry Smith       } else {
60a7e14dcfSSatish Balay         v1[i] = TAO_INFINITY;
61a7e14dcfSSatish Balay       }
62a7e14dcfSSatish Balay     }
632d0e5244SBarry Smith   } else if (2.0 == p) {
64a7e14dcfSSatish Balay     for (i = 0; i < n; ++i){
65a7e14dcfSSatish Balay       v1[i] *= v1[i];
66a7e14dcfSSatish Balay     }
672d0e5244SBarry Smith   } else if (-2.0 == p) {
68a7e14dcfSSatish Balay     for (i = 0; i < n; ++i){
69a7e14dcfSSatish Balay       v1[i] = 1.0 / (v1[i] * v1[i]);
70a7e14dcfSSatish Balay     }
712d0e5244SBarry Smith   } else {
72a7e14dcfSSatish Balay     for (i = 0; i < n; ++i) {
73a7e14dcfSSatish Balay       if (v1[i] >= 0) {
74a7e14dcfSSatish Balay         v1[i] = PetscPowScalar(v1[i], p);
752d0e5244SBarry Smith       } else {
76a7e14dcfSSatish Balay         v1[i] = TAO_INFINITY;
77a7e14dcfSSatish Balay       }
78a7e14dcfSSatish Balay     }
79a7e14dcfSSatish Balay   }
80a7e14dcfSSatish Balay   ierr = VecRestoreArray(v,&v1);CHKERRQ(ierr);
81a7e14dcfSSatish Balay   PetscFunctionReturn(0);
82a7e14dcfSSatish Balay }
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay #undef __FUNCT__
85a7e14dcfSSatish Balay #define __FUNCT__ "VecMedian"
86a7e14dcfSSatish Balay /*@
87a7e14dcfSSatish Balay   VecMedian - Computes the componentwise median of three vectors
88a7e14dcfSSatish Balay   and stores the result in this vector.  Used primarily for projecting
89a7e14dcfSSatish Balay   a vector within upper and lower bounds.
90a7e14dcfSSatish Balay 
91a7e14dcfSSatish Balay   Logically Collective
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay   Input Parameters:
94a7e14dcfSSatish Balay . Vec1, Vec2, Vec3 - The three vectors
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay   Output Parameter:
97a7e14dcfSSatish Balay . VMedian - The median vector
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay   Level: advanced
100a7e14dcfSSatish Balay @*/
101a7e14dcfSSatish Balay PetscErrorCode VecMedian(Vec Vec1, Vec Vec2, Vec Vec3, Vec VMedian)
102a7e14dcfSSatish Balay {
103a7e14dcfSSatish Balay   PetscErrorCode ierr;
104a7e14dcfSSatish Balay   PetscInt       i,n,low1,low2,low3,low4,high1,high2,high3,high4;
105a7e14dcfSSatish Balay   PetscReal      *v1,*v2,*v3,*vmed;
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay   PetscFunctionBegin;
108a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1);
109a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2);
110a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec3,VEC_CLASSID,3);
111a7e14dcfSSatish Balay   PetscValidHeaderSpecific(VMedian,VEC_CLASSID,4);
112a7e14dcfSSatish Balay 
113a7e14dcfSSatish Balay   if (Vec1==Vec2 || Vec1==Vec3){
114a7e14dcfSSatish Balay     ierr=VecCopy(Vec1,VMedian);CHKERRQ(ierr);
115a7e14dcfSSatish Balay     PetscFunctionReturn(0);
116a7e14dcfSSatish Balay   }
117a7e14dcfSSatish Balay   if (Vec2==Vec3){
118a7e14dcfSSatish Balay     ierr=VecCopy(Vec2,VMedian);CHKERRQ(ierr);
119a7e14dcfSSatish Balay     PetscFunctionReturn(0);
120a7e14dcfSSatish Balay   }
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay   PetscValidType(Vec1,1);
123a7e14dcfSSatish Balay   PetscValidType(Vec2,2);
124a7e14dcfSSatish Balay   PetscValidType(VMedian,4);
125a7e14dcfSSatish Balay   PetscCheckSameType(Vec1,1,Vec2,2); PetscCheckSameType(Vec1,1,VMedian,4);
126a7e14dcfSSatish Balay   PetscCheckSameComm(Vec1,1,Vec2,2); PetscCheckSameComm(Vec1,1,VMedian,4);
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec1, &low1, &high1);CHKERRQ(ierr);
129a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec2, &low2, &high2);CHKERRQ(ierr);
130a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec3, &low3, &high3);CHKERRQ(ierr);
131a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(VMedian, &low4, &high4);CHKERRQ(ierr);
1322d0e5244SBarry Smith   if ( low1!= low2 || low1!= low3 || low1!= low4 || high1!= high2 || high1!= high3 || high1!= high4) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths");
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay   ierr = VecGetArray(Vec1,&v1);CHKERRQ(ierr);
135a7e14dcfSSatish Balay   ierr = VecGetArray(Vec2,&v2);CHKERRQ(ierr);
136a7e14dcfSSatish Balay   ierr = VecGetArray(Vec3,&v3);CHKERRQ(ierr);
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay   if ( VMedian != Vec1 && VMedian != Vec2 && VMedian != Vec3){
139a7e14dcfSSatish Balay     ierr = VecGetArray(VMedian,&vmed);CHKERRQ(ierr);
140a7e14dcfSSatish Balay   } else if ( VMedian==Vec1 ){
141a7e14dcfSSatish Balay     vmed=v1;
142a7e14dcfSSatish Balay   } else if ( VMedian==Vec2 ){
143a7e14dcfSSatish Balay     vmed=v2;
144a7e14dcfSSatish Balay   } else {
145a7e14dcfSSatish Balay     vmed=v3;
146a7e14dcfSSatish Balay   }
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay   ierr=VecGetLocalSize(Vec1,&n);CHKERRQ(ierr);
149a7e14dcfSSatish Balay 
150a7e14dcfSSatish Balay   for (i=0;i<n;i++){
151a7e14dcfSSatish Balay     vmed[i]=PetscMax(PetscMax(PetscMin(v1[i],v2[i]),PetscMin(v1[i],v3[i])),PetscMin(v2[i],v3[i]));
152a7e14dcfSSatish Balay   }
153a7e14dcfSSatish Balay 
154a7e14dcfSSatish Balay   ierr = VecRestoreArray(Vec1,&v1);CHKERRQ(ierr);
155a7e14dcfSSatish Balay   ierr = VecRestoreArray(Vec2,&v2);CHKERRQ(ierr);
156a7e14dcfSSatish Balay   ierr = VecRestoreArray(Vec3,&v2);CHKERRQ(ierr);
157a7e14dcfSSatish Balay 
158a7e14dcfSSatish Balay   if (VMedian!=Vec1 && VMedian != Vec2 && VMedian != Vec3){
159a7e14dcfSSatish Balay     ierr = VecRestoreArray(VMedian,&vmed);CHKERRQ(ierr);
160a7e14dcfSSatish Balay   }
161a7e14dcfSSatish Balay   PetscFunctionReturn(0);
162a7e14dcfSSatish Balay }
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal Fischer(PetscReal a, PetscReal b)
165a7e14dcfSSatish Balay {
166a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
167a7e14dcfSSatish Balay    if (a + b <= 0) {
168a7e14dcfSSatish Balay      return PetscSqrtScalar(a*a + b*b) - (a + b);
169a7e14dcfSSatish Balay    }
170a7e14dcfSSatish Balay    return -2.0*a*b / (PetscSqrtScalar(a*a + b*b) + (a + b));
171a7e14dcfSSatish Balay }
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay #undef __FUNCT__
174a7e14dcfSSatish Balay #define __FUNCT__ "VecFischer"
175a7e14dcfSSatish Balay /*@
176a7e14dcfSSatish Balay    VecFischer - Evaluates the Fischer-Burmeister function for complementarity
177a7e14dcfSSatish Balay    problems.
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay    Logically Collective on vectors
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay    Input Parameters:
182a7e14dcfSSatish Balay +  X - current point
183a7e14dcfSSatish Balay .  F - function evaluated at x
184a7e14dcfSSatish Balay .  L - lower bounds
185a7e14dcfSSatish Balay -  U - upper bounds
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay    Output Parameters:
188a7e14dcfSSatish Balay .  FB - The Fischer-Burmeister function vector
189a7e14dcfSSatish Balay 
190a7e14dcfSSatish Balay    Notes:
191a7e14dcfSSatish Balay    The Fischer-Burmeister function is defined as
192a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b) - a - b
193a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
194a7e14dcfSSatish Balay    system of equations.
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay    The result of this function is done by cases:
197a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i]
198a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i])
199a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i])
200a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u]))
201a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay    Level: developer
204a7e14dcfSSatish Balay 
205a7e14dcfSSatish Balay @*/
206a7e14dcfSSatish Balay PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB)
207a7e14dcfSSatish Balay {
208a7e14dcfSSatish Balay   PetscReal      *x, *f, *l, *u, *fb;
209a7e14dcfSSatish Balay   PetscReal      xval, fval, lval, uval;
210a7e14dcfSSatish Balay   PetscErrorCode ierr;
211a7e14dcfSSatish Balay   PetscInt       low[5], high[5], n, i;
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay   PetscFunctionBegin;
214a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,1);
215a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID,2);
216a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID,3);
217a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID,4);
218a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID,4);
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr);
221a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr);
222a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr);
223a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr);
224a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr);
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
2272d0e5244SBarry Smith     if (low[0] != low[i] || high[0] != high[i]) SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
228a7e14dcfSSatish Balay   }
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay   ierr = VecGetArray(X, &x);CHKERRQ(ierr);
231a7e14dcfSSatish Balay   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
232a7e14dcfSSatish Balay   ierr = VecGetArray(L, &l);CHKERRQ(ierr);
233a7e14dcfSSatish Balay   ierr = VecGetArray(U, &u);CHKERRQ(ierr);
234a7e14dcfSSatish Balay   ierr = VecGetArray(FB, &fb);CHKERRQ(ierr);
235a7e14dcfSSatish Balay 
236a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
239a7e14dcfSSatish Balay     xval = x[i]; fval = f[i];
240a7e14dcfSSatish Balay     lval = l[i]; uval = u[i];
241a7e14dcfSSatish Balay 
242a7e14dcfSSatish Balay     if ((lval <= -TAO_INFINITY) && (uval >= TAO_INFINITY)) {
243a7e14dcfSSatish Balay       fb[i] = -fval;
2442d0e5244SBarry Smith     } else if (lval <= -TAO_INFINITY) {
245a7e14dcfSSatish Balay       fb[i] = -Fischer(uval - xval, -fval);
2462d0e5244SBarry Smith     } else if (uval >=  TAO_INFINITY) {
247a7e14dcfSSatish Balay       fb[i] =  Fischer(xval - lval,  fval);
2482d0e5244SBarry Smith     } else if (lval == uval) {
249a7e14dcfSSatish Balay       fb[i] = lval - xval;
2502d0e5244SBarry Smith     } else {
251a7e14dcfSSatish Balay       fval  =  Fischer(uval - xval, -fval);
252a7e14dcfSSatish Balay       fb[i] =  Fischer(xval - lval,  fval);
253a7e14dcfSSatish Balay     }
254a7e14dcfSSatish Balay   }
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay   ierr = VecRestoreArray(X, &x);CHKERRQ(ierr);
257a7e14dcfSSatish Balay   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
258a7e14dcfSSatish Balay   ierr = VecRestoreArray(L, &l);CHKERRQ(ierr);
259a7e14dcfSSatish Balay   ierr = VecRestoreArray(U, &u);CHKERRQ(ierr);
260a7e14dcfSSatish Balay   ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr);
261a7e14dcfSSatish Balay   PetscFunctionReturn(0);
262a7e14dcfSSatish Balay }
263a7e14dcfSSatish Balay 
264a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c)
265a7e14dcfSSatish Balay {
266a7e14dcfSSatish Balay   /* Method suggested by Bob Vanderbei */
267a7e14dcfSSatish Balay    if (a + b <= 0) {
268a7e14dcfSSatish Balay      return PetscSqrtScalar(a*a + b*b + 2.0*c*c) - (a + b);
269a7e14dcfSSatish Balay    }
270a7e14dcfSSatish Balay    return 2.0*(c*c - a*b) / (PetscSqrtScalar(a*a + b*b + 2.0*c*c) + (a + b));
271a7e14dcfSSatish Balay }
272a7e14dcfSSatish Balay 
273a7e14dcfSSatish Balay #undef __FUNCT__
274a7e14dcfSSatish Balay #define __FUNCT__ "VecSFischer"
275a7e14dcfSSatish Balay /*@
276a7e14dcfSSatish Balay    VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for
277a7e14dcfSSatish Balay    complementarity problems.
278a7e14dcfSSatish Balay 
279a7e14dcfSSatish Balay    Logically Collective on vectors
280a7e14dcfSSatish Balay 
281a7e14dcfSSatish Balay    Input Parameters:
282a7e14dcfSSatish Balay +  X - current point
283a7e14dcfSSatish Balay .  F - function evaluated at x
284a7e14dcfSSatish Balay .  L - lower bounds
285a7e14dcfSSatish Balay .  U - upper bounds
286a7e14dcfSSatish Balay -  mu - smoothing parameter
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay    Output Parameters:
289a7e14dcfSSatish Balay .  FB - The Smoothed Fischer-Burmeister function vector
290a7e14dcfSSatish Balay 
291a7e14dcfSSatish Balay    Notes:
292a7e14dcfSSatish Balay    The Smoothed Fischer-Burmeister function is defined as
293a7e14dcfSSatish Balay $        phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b
294a7e14dcfSSatish Balay    and is used reformulate a complementarity problem as a semismooth
295a7e14dcfSSatish Balay    system of equations.
296a7e14dcfSSatish Balay 
297a7e14dcfSSatish Balay    The result of this function is done by cases:
298a7e14dcfSSatish Balay +  l[i] == -infinity, u[i] == infinity  -- fb[i] = -f[i] - 2*mu*x[i]
299a7e14dcfSSatish Balay .  l[i] == -infinity, u[i] finite       -- fb[i] = phi(u[i]-x[i], -f[i], mu)
300a7e14dcfSSatish Balay .  l[i] finite,       u[i] == infinity  -- fb[i] = phi(x[i]-l[i],  f[i], mu)
301a7e14dcfSSatish Balay .  l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu)
302a7e14dcfSSatish Balay -  otherwise l[i] == u[i] -- fb[i] = l[i] - x[i]
303a7e14dcfSSatish Balay 
304a7e14dcfSSatish Balay    Level: developer
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay .seealso  VecFischer()
307a7e14dcfSSatish Balay @*/
308a7e14dcfSSatish Balay PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB)
309a7e14dcfSSatish Balay {
310a7e14dcfSSatish Balay   PetscReal      *x, *f, *l, *u, *fb;
311a7e14dcfSSatish Balay   PetscReal      xval, fval, lval, uval;
312a7e14dcfSSatish Balay   PetscErrorCode ierr;
313a7e14dcfSSatish Balay   PetscInt       low[5], high[5], n, i;
314a7e14dcfSSatish Balay 
315a7e14dcfSSatish Balay   PetscFunctionBegin;
316a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,1);
317a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F, VEC_CLASSID,2);
318a7e14dcfSSatish Balay   PetscValidHeaderSpecific(L, VEC_CLASSID,3);
319a7e14dcfSSatish Balay   PetscValidHeaderSpecific(U, VEC_CLASSID,4);
320a7e14dcfSSatish Balay   PetscValidHeaderSpecific(FB, VEC_CLASSID,6);
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr);
323a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr);
324a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr);
325a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr);
326a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr);
327a7e14dcfSSatish Balay 
328a7e14dcfSSatish Balay   for (i = 1; i < 4; ++i) {
3292d0e5244SBarry Smith     if (low[0] != low[i] || high[0] != high[i]) SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
330a7e14dcfSSatish Balay   }
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay   ierr = VecGetArray(X, &x);CHKERRQ(ierr);
333a7e14dcfSSatish Balay   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
334a7e14dcfSSatish Balay   ierr = VecGetArray(L, &l);CHKERRQ(ierr);
335a7e14dcfSSatish Balay   ierr = VecGetArray(U, &u);CHKERRQ(ierr);
336a7e14dcfSSatish Balay   ierr = VecGetArray(FB, &fb);CHKERRQ(ierr);
337a7e14dcfSSatish Balay 
338a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
339a7e14dcfSSatish Balay 
340a7e14dcfSSatish Balay   for (i = 0; i < n; ++i) {
341a7e14dcfSSatish Balay     xval = (*x++); fval = (*f++);
342a7e14dcfSSatish Balay     lval = (*l++); uval = (*u++);
343a7e14dcfSSatish Balay 
344a7e14dcfSSatish Balay     if ((lval <= -TAO_INFINITY) && (uval >= TAO_INFINITY)) {
345a7e14dcfSSatish Balay       (*fb++) = -fval - mu*xval;
3462d0e5244SBarry Smith     } else if (lval <= -TAO_INFINITY) {
347a7e14dcfSSatish Balay       (*fb++) = -SFischer(uval - xval, -fval, mu);
3482d0e5244SBarry Smith     } else if (uval >=  TAO_INFINITY) {
349a7e14dcfSSatish Balay       (*fb++) =  SFischer(xval - lval,  fval, mu);
3502d0e5244SBarry Smith     } else if (lval == uval) {
351a7e14dcfSSatish Balay       (*fb++) = lval - xval;
3522d0e5244SBarry Smith     } else {
353a7e14dcfSSatish Balay       fval    =  SFischer(uval - xval, -fval, mu);
354a7e14dcfSSatish Balay       (*fb++) =  SFischer(xval - lval,  fval, mu);
355a7e14dcfSSatish Balay     }
356a7e14dcfSSatish Balay   }
357a7e14dcfSSatish Balay   x -= n; f -= n; l -=n; u -= n; fb -= n;
358a7e14dcfSSatish Balay 
359a7e14dcfSSatish Balay   ierr = VecRestoreArray(X, &x);CHKERRQ(ierr);
360a7e14dcfSSatish Balay   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
361a7e14dcfSSatish Balay   ierr = VecRestoreArray(L, &l);CHKERRQ(ierr);
362a7e14dcfSSatish Balay   ierr = VecRestoreArray(U, &u);CHKERRQ(ierr);
363a7e14dcfSSatish Balay   ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr);
364a7e14dcfSSatish Balay   PetscFunctionReturn(0);
365a7e14dcfSSatish Balay }
366a7e14dcfSSatish Balay 
367a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal fischnorm(PetscReal a, PetscReal b)
368a7e14dcfSSatish Balay {
369a7e14dcfSSatish Balay   return PetscSqrtScalar(a*a + b*b);
370a7e14dcfSSatish Balay }
371a7e14dcfSSatish Balay 
372a7e14dcfSSatish Balay PETSC_STATIC_INLINE PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c)
373a7e14dcfSSatish Balay {
374a7e14dcfSSatish Balay   return PetscSqrtScalar(a*a + b*b + 2.0*c*c);
375a7e14dcfSSatish Balay }
376a7e14dcfSSatish Balay 
377a7e14dcfSSatish Balay #undef __FUNCT__
378a7e14dcfSSatish Balay #define __FUNCT__ "D_Fischer"
379a7e14dcfSSatish Balay /*@
380a7e14dcfSSatish Balay    D_Fischer - Calculates an element of the B-subdifferential of the
381a7e14dcfSSatish Balay    Fischer-Burmeister function for complementarity problems.
382a7e14dcfSSatish Balay 
383a7e14dcfSSatish Balay    Collective on jac
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay    Input Parameters:
386a7e14dcfSSatish Balay +  jac - the jacobian of f at X
387a7e14dcfSSatish Balay .  X - current point
388a7e14dcfSSatish Balay .  Con - constraints function evaluated at X
389a7e14dcfSSatish Balay .  XL - lower bounds
390a7e14dcfSSatish Balay .  XU - upper bounds
391a7e14dcfSSatish Balay .  t1 - work vector
392a7e14dcfSSatish Balay -  t2 - work vector
393a7e14dcfSSatish Balay 
394a7e14dcfSSatish Balay    Output Parameters:
395a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
396a7e14dcfSSatish Balay -  Db - row scaling component of the result
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay    Level: developer
399a7e14dcfSSatish Balay 
400a7e14dcfSSatish Balay .seealso: VecFischer()
401a7e14dcfSSatish Balay @*/
4022d0e5244SBarry Smith PetscErrorCode D_Fischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db)
403a7e14dcfSSatish Balay {
404a7e14dcfSSatish Balay   PetscErrorCode ierr;
405a7e14dcfSSatish Balay   PetscInt       i,nn;
406a7e14dcfSSatish Balay   PetscReal      *x,*f,*l,*u,*da,*db,*t1,*t2;
407a7e14dcfSSatish Balay   PetscReal      ai,bi,ci,di,ei;
408a7e14dcfSSatish Balay 
409a7e14dcfSSatish Balay   PetscFunctionBegin;
410a7e14dcfSSatish Balay   ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr);
411a7e14dcfSSatish Balay   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
412a7e14dcfSSatish Balay   ierr = VecGetArray(Con,&f);CHKERRQ(ierr);
413a7e14dcfSSatish Balay   ierr = VecGetArray(XL,&l);CHKERRQ(ierr);
414a7e14dcfSSatish Balay   ierr = VecGetArray(XU,&u);CHKERRQ(ierr);
415a7e14dcfSSatish Balay   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
416a7e14dcfSSatish Balay   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
417a7e14dcfSSatish Balay   ierr = VecGetArray(T1,&t1);CHKERRQ(ierr);
418a7e14dcfSSatish Balay   ierr = VecGetArray(T2,&t2);CHKERRQ(ierr);
419a7e14dcfSSatish Balay 
420a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
421a7e14dcfSSatish Balay     da[i] = 0;
422a7e14dcfSSatish Balay     db[i] = 0;
423a7e14dcfSSatish Balay     t1[i] = 0;
424a7e14dcfSSatish Balay 
425a7e14dcfSSatish Balay     if (PetscAbsReal(f[i]) <= PETSC_MACHINE_EPSILON) {
426a7e14dcfSSatish Balay       if (l[i] > TAO_NINFINITY && PetscAbsReal(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) {
427a7e14dcfSSatish Balay         t1[i] = 1;
428a7e14dcfSSatish Balay         da[i] = 1;
429a7e14dcfSSatish Balay       }
430a7e14dcfSSatish Balay 
431a7e14dcfSSatish Balay       if (u[i] <  TAO_INFINITY && PetscAbsReal(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) {
432a7e14dcfSSatish Balay         t1[i] = 1;
433a7e14dcfSSatish Balay         db[i] = 1;
434a7e14dcfSSatish Balay       }
435a7e14dcfSSatish Balay     }
436a7e14dcfSSatish Balay   }
437a7e14dcfSSatish Balay 
438a7e14dcfSSatish Balay   ierr = VecRestoreArray(T1,&t1);CHKERRQ(ierr);
439a7e14dcfSSatish Balay   ierr = VecRestoreArray(T2,&t2);CHKERRQ(ierr);
440a7e14dcfSSatish Balay   ierr = MatMult(jac,T1,T2);CHKERRQ(ierr);
441a7e14dcfSSatish Balay   ierr = VecGetArray(T2,&t2);CHKERRQ(ierr);
442a7e14dcfSSatish Balay 
443a7e14dcfSSatish Balay   for (i = 0; i < nn; i++) {
444a7e14dcfSSatish Balay     if ((l[i] <= TAO_NINFINITY) && (u[i] >= TAO_INFINITY)) {
445a7e14dcfSSatish Balay       da[i] = 0;
446a7e14dcfSSatish Balay       db[i] = -1;
4472d0e5244SBarry Smith     } else if (l[i] <= TAO_NINFINITY) {
448a7e14dcfSSatish Balay       if (db[i] >= 1) {
449a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
450a7e14dcfSSatish Balay 
451a7e14dcfSSatish Balay         da[i] = -1/ai - 1;
452a7e14dcfSSatish Balay         db[i] = -t2[i]/ai - 1;
4532d0e5244SBarry Smith       } else {
454a7e14dcfSSatish Balay         bi = u[i] - x[i];
455a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
456a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
457a7e14dcfSSatish Balay 
458a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
459a7e14dcfSSatish Balay         db[i] = -f[i] / ai - 1;
460a7e14dcfSSatish Balay       }
4612d0e5244SBarry Smith     } else if (u[i] >=  TAO_INFINITY) {
462a7e14dcfSSatish Balay       if (da[i] >= 1) {
463a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
464a7e14dcfSSatish Balay 
465a7e14dcfSSatish Balay         da[i] = 1 / ai - 1;
466a7e14dcfSSatish Balay         db[i] = t2[i] / ai - 1;
4672d0e5244SBarry Smith       } else {
468a7e14dcfSSatish Balay         bi = x[i] - l[i];
469a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
470a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
471a7e14dcfSSatish Balay 
472a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
473a7e14dcfSSatish Balay         db[i] = f[i] / ai - 1;
474a7e14dcfSSatish Balay       }
4752d0e5244SBarry Smith     } else if (l[i] == u[i]) {
476a7e14dcfSSatish Balay       da[i] = -1;
477a7e14dcfSSatish Balay       db[i] = 0;
4782d0e5244SBarry Smith     } else {
479a7e14dcfSSatish Balay       if (db[i] >= 1) {
480a7e14dcfSSatish Balay         ai = fischnorm(1, t2[i]);
481a7e14dcfSSatish Balay 
482a7e14dcfSSatish Balay         ci = 1 / ai + 1;
483a7e14dcfSSatish Balay         di = t2[i] / ai + 1;
4842d0e5244SBarry Smith       } else {
485a7e14dcfSSatish Balay         bi = x[i] - u[i];
486a7e14dcfSSatish Balay         ai = fischnorm(bi, f[i]);
487a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
488a7e14dcfSSatish Balay 
489a7e14dcfSSatish Balay         ci = bi / ai + 1;
490a7e14dcfSSatish Balay         di = f[i] / ai + 1;
491a7e14dcfSSatish Balay       }
492a7e14dcfSSatish Balay 
493a7e14dcfSSatish Balay       if (da[i] >= 1) {
494a7e14dcfSSatish Balay         bi = ci + di*t2[i];
495a7e14dcfSSatish Balay         ai = fischnorm(1, bi);
496a7e14dcfSSatish Balay 
497a7e14dcfSSatish Balay         bi = bi / ai - 1;
498a7e14dcfSSatish Balay         ai = 1 / ai - 1;
4992d0e5244SBarry Smith       } else {
500a7e14dcfSSatish Balay         ei = Fischer(u[i] - x[i], -f[i]);
501a7e14dcfSSatish Balay         ai = fischnorm(x[i] - l[i], ei);
502a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
503a7e14dcfSSatish Balay 
504a7e14dcfSSatish Balay         bi = ei / ai - 1;
505a7e14dcfSSatish Balay         ai = (x[i] - l[i]) / ai - 1;
506a7e14dcfSSatish Balay       }
507a7e14dcfSSatish Balay 
508a7e14dcfSSatish Balay       da[i] = ai + bi*ci;
509a7e14dcfSSatish Balay       db[i] = bi*di;
510a7e14dcfSSatish Balay     }
511a7e14dcfSSatish Balay   }
512a7e14dcfSSatish Balay 
513a7e14dcfSSatish Balay   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
514a7e14dcfSSatish Balay   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
515a7e14dcfSSatish Balay   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
516a7e14dcfSSatish Balay   ierr = VecRestoreArray(Con,&f);CHKERRQ(ierr);
517a7e14dcfSSatish Balay   ierr = VecRestoreArray(XL,&l);CHKERRQ(ierr);
518a7e14dcfSSatish Balay   ierr = VecRestoreArray(XU,&u);CHKERRQ(ierr);
519a7e14dcfSSatish Balay   ierr = VecRestoreArray(T2,&t2);CHKERRQ(ierr);
520a7e14dcfSSatish Balay   PetscFunctionReturn(0);
5218e3154b5SSatish Balay }
522a7e14dcfSSatish Balay 
523a7e14dcfSSatish Balay #undef __FUNCT__
524a7e14dcfSSatish Balay #define __FUNCT__ "D_SFischer"
525a7e14dcfSSatish Balay /*@
526a7e14dcfSSatish Balay    D_SFischer - Calculates an element of the B-subdifferential of the
527a7e14dcfSSatish Balay    smoothed Fischer-Burmeister function for complementarity problems.
528a7e14dcfSSatish Balay 
529a7e14dcfSSatish Balay    Collective on jac
530a7e14dcfSSatish Balay 
531a7e14dcfSSatish Balay    Input Parameters:
532a7e14dcfSSatish Balay +  jac - the jacobian of f at X
533a7e14dcfSSatish Balay .  X - current point
534a7e14dcfSSatish Balay .  F - constraint function evaluated at X
535a7e14dcfSSatish Balay .  XL - lower bounds
536a7e14dcfSSatish Balay .  XU - upper bounds
537a7e14dcfSSatish Balay .  mu - smoothing parameter
538a7e14dcfSSatish Balay .  T1 - work vector
539a7e14dcfSSatish Balay -  T2 - work vector
540a7e14dcfSSatish Balay 
541a7e14dcfSSatish Balay    Output Parameter:
542a7e14dcfSSatish Balay +  Da - diagonal perturbation component of the result
543a7e14dcfSSatish Balay .  Db - row scaling component of the result
544a7e14dcfSSatish Balay -  Dm - derivative with respect to scaling parameter
545a7e14dcfSSatish Balay 
546a7e14dcfSSatish Balay    Level: developer
547a7e14dcfSSatish Balay 
548a7e14dcfSSatish Balay .seealso D_Fischer()
549a7e14dcfSSatish Balay @*/
5502d0e5244SBarry 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)
551a7e14dcfSSatish Balay {
552a7e14dcfSSatish Balay   PetscErrorCode ierr;
553a7e14dcfSSatish Balay   PetscInt       i,nn;
554a7e14dcfSSatish Balay   PetscReal      *x, *f, *l, *u, *da, *db, *dm;
555a7e14dcfSSatish Balay   PetscReal      ai, bi, ci, di, ei, fi;
556a7e14dcfSSatish Balay 
557a7e14dcfSSatish Balay   PetscFunctionBegin;
558a7e14dcfSSatish Balay   if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) {
559a7e14dcfSSatish Balay     ierr = VecZeroEntries(Dm);CHKERRQ(ierr);
560a7e14dcfSSatish Balay     ierr = D_Fischer(jac, X, Con, XL, XU, T1, T2, Da, Db);CHKERRQ(ierr);
5612d0e5244SBarry Smith   } else {
562a7e14dcfSSatish Balay     ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr);
563a7e14dcfSSatish Balay     ierr = VecGetArray(X,&x);CHKERRQ(ierr);
564a7e14dcfSSatish Balay     ierr = VecGetArray(Con,&f);CHKERRQ(ierr);
565a7e14dcfSSatish Balay     ierr = VecGetArray(XL,&l);CHKERRQ(ierr);
566a7e14dcfSSatish Balay     ierr = VecGetArray(XU,&u);CHKERRQ(ierr);
567a7e14dcfSSatish Balay     ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
568a7e14dcfSSatish Balay     ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
569a7e14dcfSSatish Balay     ierr = VecGetArray(Dm,&dm);CHKERRQ(ierr);
570a7e14dcfSSatish Balay 
571a7e14dcfSSatish Balay     for (i = 0; i < nn; ++i) {
572a7e14dcfSSatish Balay       if ((l[i] <= TAO_NINFINITY) && (u[i] >= TAO_INFINITY)) {
573a7e14dcfSSatish Balay         da[i] = -mu;
574a7e14dcfSSatish Balay         db[i] = -1;
575a7e14dcfSSatish Balay         dm[i] = -x[i];
5762d0e5244SBarry Smith       } else if (l[i] <= TAO_NINFINITY) {
577a7e14dcfSSatish Balay         bi = u[i] - x[i];
578a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
579a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
580a7e14dcfSSatish Balay 
581a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
582a7e14dcfSSatish Balay         db[i] = -f[i] / ai - 1;
583a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
5842d0e5244SBarry Smith       } else if (u[i] >=  TAO_INFINITY) {
585a7e14dcfSSatish Balay         bi = x[i] - l[i];
586a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
587a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
588a7e14dcfSSatish Balay 
589a7e14dcfSSatish Balay         da[i] = bi / ai - 1;
590a7e14dcfSSatish Balay         db[i] = f[i] / ai - 1;
591a7e14dcfSSatish Balay         dm[i] = 2.0 * mu / ai;
5922d0e5244SBarry Smith       } else if (l[i] == u[i]) {
593a7e14dcfSSatish Balay         da[i] = -1;
594a7e14dcfSSatish Balay         db[i] = 0;
595a7e14dcfSSatish Balay         dm[i] = 0;
5962d0e5244SBarry Smith       } else {
597a7e14dcfSSatish Balay         bi = x[i] - u[i];
598a7e14dcfSSatish Balay         ai = fischsnorm(bi, f[i], mu);
599a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay         ci = bi / ai + 1;
602a7e14dcfSSatish Balay         di = f[i] / ai + 1;
603a7e14dcfSSatish Balay         fi = 2.0 * mu / ai;
604a7e14dcfSSatish Balay 
605a7e14dcfSSatish Balay         ei = SFischer(u[i] - x[i], -f[i], mu);
606a7e14dcfSSatish Balay         ai = fischsnorm(x[i] - l[i], ei, mu);
607a7e14dcfSSatish Balay         ai = PetscMax(PETSC_MACHINE_EPSILON, ai);
608a7e14dcfSSatish Balay 
609a7e14dcfSSatish Balay         bi = ei / ai - 1;
610a7e14dcfSSatish Balay         ei = 2.0 * mu / ei;
611a7e14dcfSSatish Balay         ai = (x[i] - l[i]) / ai - 1;
612a7e14dcfSSatish Balay 
613a7e14dcfSSatish Balay         da[i] = ai + bi*ci;
614a7e14dcfSSatish Balay         db[i] = bi*di;
615a7e14dcfSSatish Balay         dm[i] = ei + bi*fi;
616a7e14dcfSSatish Balay       }
617a7e14dcfSSatish Balay     }
618a7e14dcfSSatish Balay 
619a7e14dcfSSatish Balay     ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
620a7e14dcfSSatish Balay     ierr = VecRestoreArray(Con,&f);CHKERRQ(ierr);
621a7e14dcfSSatish Balay     ierr = VecRestoreArray(XL,&l);CHKERRQ(ierr);
622a7e14dcfSSatish Balay     ierr = VecRestoreArray(XU,&u);CHKERRQ(ierr);
623a7e14dcfSSatish Balay     ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
624a7e14dcfSSatish Balay     ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
625a7e14dcfSSatish Balay     ierr = VecRestoreArray(Dm,&dm);CHKERRQ(ierr);
626a7e14dcfSSatish Balay   }
627a7e14dcfSSatish Balay   PetscFunctionReturn(0);
628a7e14dcfSSatish Balay }
629a7e14dcfSSatish Balay 
630