1af0996ceSBarry Smith #include <petsc/private/petscimpl.h> 2ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 38370d7cdSHansol Suh #include <petscsys.h> 4a7e14dcfSSatish Balay 5*9371c9d4SSatish Balay static inline PetscReal Fischer(PetscReal a, PetscReal b) { 6a7e14dcfSSatish Balay /* Method suggested by Bob Vanderbei */ 7*9371c9d4SSatish Balay if (a + b <= 0) { return PetscSqrtReal(a * a + b * b) - (a + b); } 846bdf8c8SLisandro Dalcin return -2.0 * a * b / (PetscSqrtReal(a * a + b * b) + (a + b)); 9a7e14dcfSSatish Balay } 10a7e14dcfSSatish Balay 11a7e14dcfSSatish Balay /*@ 12a7e14dcfSSatish Balay VecFischer - Evaluates the Fischer-Burmeister function for complementarity 13a7e14dcfSSatish Balay problems. 14a7e14dcfSSatish Balay 15a7e14dcfSSatish Balay Logically Collective on vectors 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay Input Parameters: 18a7e14dcfSSatish Balay + X - current point 19a7e14dcfSSatish Balay . F - function evaluated at x 20a7e14dcfSSatish Balay . L - lower bounds 21a7e14dcfSSatish Balay - U - upper bounds 22a7e14dcfSSatish Balay 23f899ff85SJose E. Roman Output Parameter: 24a7e14dcfSSatish Balay . FB - The Fischer-Burmeister function vector 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay Notes: 27a7e14dcfSSatish Balay The Fischer-Burmeister function is defined as 28a7e14dcfSSatish Balay $ phi(a,b) := sqrt(a*a + b*b) - a - b 29a7e14dcfSSatish Balay and is used reformulate a complementarity problem as a semismooth 30a7e14dcfSSatish Balay system of equations. 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay The result of this function is done by cases: 33a7e14dcfSSatish Balay + l[i] == -infinity, u[i] == infinity -- fb[i] = -f[i] 34a7e14dcfSSatish Balay . l[i] == -infinity, u[i] finite -- fb[i] = phi(u[i]-x[i], -f[i]) 35a7e14dcfSSatish Balay . l[i] finite, u[i] == infinity -- fb[i] = phi(x[i]-l[i], f[i]) 36a7e14dcfSSatish Balay . l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u])) 37a7e14dcfSSatish Balay - otherwise l[i] == u[i] -- fb[i] = l[i] - x[i] 38a7e14dcfSSatish Balay 39a7e14dcfSSatish Balay Level: developer 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay @*/ 42*9371c9d4SSatish Balay PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB) { 4346bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 4446bdf8c8SLisandro Dalcin PetscScalar *fb; 45a7e14dcfSSatish Balay PetscReal xval, fval, lval, uval; 46a7e14dcfSSatish Balay PetscInt low[5], high[5], n, i; 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay PetscFunctionBegin; 49a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 50a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID, 2); 5176be6f4fSStefano Zampini if (L) PetscValidHeaderSpecific(L, VEC_CLASSID, 3); 5276be6f4fSStefano Zampini if (U) PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 53064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(FB, VEC_CLASSID, 5); 54a7e14dcfSSatish Balay 5576be6f4fSStefano Zampini if (!L && !U) { 5676be6f4fSStefano Zampini PetscCall(VecAXPBY(FB, -1.0, 0.0, F)); 5776be6f4fSStefano Zampini PetscFunctionReturn(0); 5876be6f4fSStefano Zampini } 5976be6f4fSStefano Zampini 609566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, low, high)); 619566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(F, low + 1, high + 1)); 629566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(L, low + 2, high + 2)); 639566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(U, low + 3, high + 3)); 649566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4)); 65a7e14dcfSSatish Balay 66*9371c9d4SSatish Balay for (i = 1; i < 4; ++i) { PetscCheck(low[0] == low[i] && high[0] == high[i], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Vectors must be identically loaded over processors"); } 67a7e14dcfSSatish Balay 689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(L, &l)); 719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArray(FB, &fb)); 73a7e14dcfSSatish Balay 749566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay for (i = 0; i < n; ++i) { 7776be6f4fSStefano Zampini xval = PetscRealPart(x[i]); 7876be6f4fSStefano Zampini fval = PetscRealPart(f[i]); 7976be6f4fSStefano Zampini lval = PetscRealPart(l[i]); 8076be6f4fSStefano Zampini uval = PetscRealPart(u[i]); 81a7e14dcfSSatish Balay 8276be6f4fSStefano Zampini if (lval <= -PETSC_INFINITY && uval >= PETSC_INFINITY) { 83a7e14dcfSSatish Balay fb[i] = -fval; 84e270355aSBarry Smith } else if (lval <= -PETSC_INFINITY) { 85a7e14dcfSSatish Balay fb[i] = -Fischer(uval - xval, -fval); 86e270355aSBarry Smith } else if (uval >= PETSC_INFINITY) { 87a7e14dcfSSatish Balay fb[i] = Fischer(xval - lval, fval); 882d0e5244SBarry Smith } else if (lval == uval) { 89a7e14dcfSSatish Balay fb[i] = lval - xval; 902d0e5244SBarry Smith } else { 91a7e14dcfSSatish Balay fval = Fischer(uval - xval, -fval); 92a7e14dcfSSatish Balay fb[i] = Fischer(xval - lval, fval); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay } 95a7e14dcfSSatish Balay 969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(L, &l)); 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(FB, &fb)); 101a7e14dcfSSatish Balay PetscFunctionReturn(0); 102a7e14dcfSSatish Balay } 103a7e14dcfSSatish Balay 104*9371c9d4SSatish Balay static inline PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c) { 105a7e14dcfSSatish Balay /* Method suggested by Bob Vanderbei */ 106*9371c9d4SSatish Balay if (a + b <= 0) { return PetscSqrtReal(a * a + b * b + 2.0 * c * c) - (a + b); } 1073f6ba705SLisandro Dalcin return 2.0 * (c * c - a * b) / (PetscSqrtReal(a * a + b * b + 2.0 * c * c) + (a + b)); 108a7e14dcfSSatish Balay } 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay /*@ 111a7e14dcfSSatish Balay VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for 112a7e14dcfSSatish Balay complementarity problems. 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay Logically Collective on vectors 115a7e14dcfSSatish Balay 116a7e14dcfSSatish Balay Input Parameters: 117a7e14dcfSSatish Balay + X - current point 118a7e14dcfSSatish Balay . F - function evaluated at x 119a7e14dcfSSatish Balay . L - lower bounds 120a7e14dcfSSatish Balay . U - upper bounds 121a7e14dcfSSatish Balay - mu - smoothing parameter 122a7e14dcfSSatish Balay 123f899ff85SJose E. Roman Output Parameter: 124a7e14dcfSSatish Balay . FB - The Smoothed Fischer-Burmeister function vector 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay Notes: 127a7e14dcfSSatish Balay The Smoothed Fischer-Burmeister function is defined as 128a7e14dcfSSatish Balay $ phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b 129a7e14dcfSSatish Balay and is used reformulate a complementarity problem as a semismooth 130a7e14dcfSSatish Balay system of equations. 131a7e14dcfSSatish Balay 132a7e14dcfSSatish Balay The result of this function is done by cases: 133a7e14dcfSSatish Balay + l[i] == -infinity, u[i] == infinity -- fb[i] = -f[i] - 2*mu*x[i] 134a7e14dcfSSatish Balay . l[i] == -infinity, u[i] finite -- fb[i] = phi(u[i]-x[i], -f[i], mu) 135a7e14dcfSSatish Balay . l[i] finite, u[i] == infinity -- fb[i] = phi(x[i]-l[i], f[i], mu) 136a7e14dcfSSatish Balay . l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu) 137a7e14dcfSSatish Balay - otherwise l[i] == u[i] -- fb[i] = l[i] - x[i] 138a7e14dcfSSatish Balay 139a7e14dcfSSatish Balay Level: developer 140a7e14dcfSSatish Balay 141db781477SPatrick Sanan .seealso `VecFischer()` 142a7e14dcfSSatish Balay @*/ 143*9371c9d4SSatish Balay PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB) { 14446bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 14546bdf8c8SLisandro Dalcin PetscScalar *fb; 146a7e14dcfSSatish Balay PetscReal xval, fval, lval, uval; 147a7e14dcfSSatish Balay PetscInt low[5], high[5], n, i; 148a7e14dcfSSatish Balay 149a7e14dcfSSatish Balay PetscFunctionBegin; 150a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 151a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID, 2); 152a7e14dcfSSatish Balay PetscValidHeaderSpecific(L, VEC_CLASSID, 3); 153a7e14dcfSSatish Balay PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 154a7e14dcfSSatish Balay PetscValidHeaderSpecific(FB, VEC_CLASSID, 6); 155a7e14dcfSSatish Balay 1569566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, low, high)); 1579566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(F, low + 1, high + 1)); 1589566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(L, low + 2, high + 2)); 1599566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(U, low + 3, high + 3)); 1609566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4)); 161a7e14dcfSSatish Balay 162*9371c9d4SSatish Balay for (i = 1; i < 4; ++i) { PetscCheck(low[0] == low[i] && high[0] == high[i], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Vectors must be identically loaded over processors"); } 163a7e14dcfSSatish Balay 1649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(L, &l)); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1689566063dSJacob Faibussowitsch PetscCall(VecGetArray(FB, &fb)); 169a7e14dcfSSatish Balay 1709566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 171a7e14dcfSSatish Balay 172a7e14dcfSSatish Balay for (i = 0; i < n; ++i) { 173*9371c9d4SSatish Balay xval = PetscRealPart(*x++); 174*9371c9d4SSatish Balay fval = PetscRealPart(*f++); 175*9371c9d4SSatish Balay lval = PetscRealPart(*l++); 176*9371c9d4SSatish Balay uval = PetscRealPart(*u++); 177a7e14dcfSSatish Balay 178e270355aSBarry Smith if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) { 179a7e14dcfSSatish Balay (*fb++) = -fval - mu * xval; 180e270355aSBarry Smith } else if (lval <= -PETSC_INFINITY) { 181a7e14dcfSSatish Balay (*fb++) = -SFischer(uval - xval, -fval, mu); 182e270355aSBarry Smith } else if (uval >= PETSC_INFINITY) { 183a7e14dcfSSatish Balay (*fb++) = SFischer(xval - lval, fval, mu); 1842d0e5244SBarry Smith } else if (lval == uval) { 185a7e14dcfSSatish Balay (*fb++) = lval - xval; 1862d0e5244SBarry Smith } else { 187a7e14dcfSSatish Balay fval = SFischer(uval - xval, -fval, mu); 188a7e14dcfSSatish Balay (*fb++) = SFischer(xval - lval, fval, mu); 189a7e14dcfSSatish Balay } 190a7e14dcfSSatish Balay } 191*9371c9d4SSatish Balay x -= n; 192*9371c9d4SSatish Balay f -= n; 193*9371c9d4SSatish Balay l -= n; 194*9371c9d4SSatish Balay u -= n; 195*9371c9d4SSatish Balay fb -= n; 196a7e14dcfSSatish Balay 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 1999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(L, &l)); 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(FB, &fb)); 202a7e14dcfSSatish Balay PetscFunctionReturn(0); 203a7e14dcfSSatish Balay } 204a7e14dcfSSatish Balay 205*9371c9d4SSatish Balay static inline PetscReal fischnorm(PetscReal a, PetscReal b) { 206658c1fc4SLisandro Dalcin return PetscSqrtReal(a * a + b * b); 207a7e14dcfSSatish Balay } 208a7e14dcfSSatish Balay 209*9371c9d4SSatish Balay static inline PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c) { 210658c1fc4SLisandro Dalcin return PetscSqrtReal(a * a + b * b + 2.0 * c * c); 211a7e14dcfSSatish Balay } 212a7e14dcfSSatish Balay 213a7e14dcfSSatish Balay /*@ 214235fd6e6SBarry Smith MatDFischer - Calculates an element of the B-subdifferential of the 215a7e14dcfSSatish Balay Fischer-Burmeister function for complementarity problems. 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay Collective on jac 218a7e14dcfSSatish Balay 219a7e14dcfSSatish Balay Input Parameters: 220a7e14dcfSSatish Balay + jac - the jacobian of f at X 221a7e14dcfSSatish Balay . X - current point 222a7e14dcfSSatish Balay . Con - constraints function evaluated at X 223a7e14dcfSSatish Balay . XL - lower bounds 224a7e14dcfSSatish Balay . XU - upper bounds 225a7e14dcfSSatish Balay . t1 - work vector 226a7e14dcfSSatish Balay - t2 - work vector 227a7e14dcfSSatish Balay 228a7e14dcfSSatish Balay Output Parameters: 229a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result 230a7e14dcfSSatish Balay - Db - row scaling component of the result 231a7e14dcfSSatish Balay 232a7e14dcfSSatish Balay Level: developer 233a7e14dcfSSatish Balay 234db781477SPatrick Sanan .seealso: `VecFischer()` 235a7e14dcfSSatish Balay @*/ 236*9371c9d4SSatish Balay PetscErrorCode MatDFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db) { 237a7e14dcfSSatish Balay PetscInt i, nn; 23846bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u, *t2; 23946bdf8c8SLisandro Dalcin PetscScalar *da, *db, *t1; 240a7e14dcfSSatish Balay PetscReal ai, bi, ci, di, ei; 241a7e14dcfSSatish Balay 242a7e14dcfSSatish Balay PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nn)); 2449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Con, &f)); 2469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &l)); 2479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &u)); 2489566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 2499566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 2509566063dSJacob Faibussowitsch PetscCall(VecGetArray(T1, &t1)); 2519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(T2, &t2)); 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay for (i = 0; i < nn; i++) { 25446bdf8c8SLisandro Dalcin da[i] = 0.0; 25546bdf8c8SLisandro Dalcin db[i] = 0.0; 25646bdf8c8SLisandro Dalcin t1[i] = 0.0; 257a7e14dcfSSatish Balay 25846bdf8c8SLisandro Dalcin if (PetscAbsScalar(f[i]) <= PETSC_MACHINE_EPSILON) { 25946bdf8c8SLisandro Dalcin if (PetscRealPart(l[i]) > PETSC_NINFINITY && PetscAbsScalar(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) { 26046bdf8c8SLisandro Dalcin t1[i] = 1.0; 26146bdf8c8SLisandro Dalcin da[i] = 1.0; 262a7e14dcfSSatish Balay } 263a7e14dcfSSatish Balay 26446bdf8c8SLisandro Dalcin if (PetscRealPart(u[i]) < PETSC_INFINITY && PetscAbsScalar(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) { 26546bdf8c8SLisandro Dalcin t1[i] = 1.0; 26646bdf8c8SLisandro Dalcin db[i] = 1.0; 267a7e14dcfSSatish Balay } 268a7e14dcfSSatish Balay } 269a7e14dcfSSatish Balay } 270a7e14dcfSSatish Balay 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(T1, &t1)); 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(T2, &t2)); 2739566063dSJacob Faibussowitsch PetscCall(MatMult(jac, T1, T2)); 2749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(T2, &t2)); 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay for (i = 0; i < nn; i++) { 27746bdf8c8SLisandro Dalcin if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { 27846bdf8c8SLisandro Dalcin da[i] = 0.0; 27946bdf8c8SLisandro Dalcin db[i] = -1.0; 28046bdf8c8SLisandro Dalcin } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { 28146bdf8c8SLisandro Dalcin if (PetscRealPart(db[i]) >= 1) { 282658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 283a7e14dcfSSatish Balay 28446bdf8c8SLisandro Dalcin da[i] = -1.0 / ai - 1.0; 28546bdf8c8SLisandro Dalcin db[i] = -t2[i] / ai - 1.0; 2862d0e5244SBarry Smith } else { 287658c1fc4SLisandro Dalcin bi = PetscRealPart(u[i]) - PetscRealPart(x[i]); 288658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 289a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 290a7e14dcfSSatish Balay 29146bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 29246bdf8c8SLisandro Dalcin db[i] = -f[i] / ai - 1.0; 293a7e14dcfSSatish Balay } 29446bdf8c8SLisandro Dalcin } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { 29546bdf8c8SLisandro Dalcin if (PetscRealPart(da[i]) >= 1) { 296658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 297a7e14dcfSSatish Balay 29846bdf8c8SLisandro Dalcin da[i] = 1.0 / ai - 1.0; 29946bdf8c8SLisandro Dalcin db[i] = t2[i] / ai - 1.0; 3002d0e5244SBarry Smith } else { 301658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(l[i]); 302658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 303a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 304a7e14dcfSSatish Balay 30546bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 30646bdf8c8SLisandro Dalcin db[i] = f[i] / ai - 1.0; 307a7e14dcfSSatish Balay } 308658c1fc4SLisandro Dalcin } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) { 30946bdf8c8SLisandro Dalcin da[i] = -1.0; 31046bdf8c8SLisandro Dalcin db[i] = 0.0; 3112d0e5244SBarry Smith } else { 31246bdf8c8SLisandro Dalcin if (PetscRealPart(db[i]) >= 1) { 313658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 314a7e14dcfSSatish Balay 31546bdf8c8SLisandro Dalcin ci = 1.0 / ai + 1.0; 316658c1fc4SLisandro Dalcin di = PetscRealPart(t2[i]) / ai + 1.0; 3172d0e5244SBarry Smith } else { 318658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(u[i]); 319658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 320a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 321a7e14dcfSSatish Balay 32246bdf8c8SLisandro Dalcin ci = bi / ai + 1.0; 323658c1fc4SLisandro Dalcin di = PetscRealPart(f[i]) / ai + 1.0; 324a7e14dcfSSatish Balay } 325a7e14dcfSSatish Balay 32646bdf8c8SLisandro Dalcin if (PetscRealPart(da[i]) >= 1) { 327658c1fc4SLisandro Dalcin bi = ci + di * PetscRealPart(t2[i]); 328658c1fc4SLisandro Dalcin ai = fischnorm(1.0, bi); 329a7e14dcfSSatish Balay 33046bdf8c8SLisandro Dalcin bi = bi / ai - 1.0; 33146bdf8c8SLisandro Dalcin ai = 1.0 / ai - 1.0; 3322d0e5244SBarry Smith } else { 333658c1fc4SLisandro Dalcin ei = Fischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i])); 334658c1fc4SLisandro Dalcin ai = fischnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei); 335a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 336a7e14dcfSSatish Balay 33746bdf8c8SLisandro Dalcin bi = ei / ai - 1.0; 338658c1fc4SLisandro Dalcin ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0; 339a7e14dcfSSatish Balay } 340a7e14dcfSSatish Balay 341a7e14dcfSSatish Balay da[i] = ai + bi * ci; 342a7e14dcfSSatish Balay db[i] = bi * di; 343a7e14dcfSSatish Balay } 344a7e14dcfSSatish Balay } 345a7e14dcfSSatish Balay 3469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 3479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 3489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Con, &f)); 3509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &l)); 3519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &u)); 3529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(T2, &t2)); 353a7e14dcfSSatish Balay PetscFunctionReturn(0); 3548e3154b5SSatish Balay } 355a7e14dcfSSatish Balay 356a7e14dcfSSatish Balay /*@ 357235fd6e6SBarry Smith MatDSFischer - Calculates an element of the B-subdifferential of the 358a7e14dcfSSatish Balay smoothed Fischer-Burmeister function for complementarity problems. 359a7e14dcfSSatish Balay 360a7e14dcfSSatish Balay Collective on jac 361a7e14dcfSSatish Balay 362a7e14dcfSSatish Balay Input Parameters: 363a7e14dcfSSatish Balay + jac - the jacobian of f at X 364a7e14dcfSSatish Balay . X - current point 365a7e14dcfSSatish Balay . F - constraint function evaluated at X 366a7e14dcfSSatish Balay . XL - lower bounds 367a7e14dcfSSatish Balay . XU - upper bounds 368a7e14dcfSSatish Balay . mu - smoothing parameter 369a7e14dcfSSatish Balay . T1 - work vector 370a7e14dcfSSatish Balay - T2 - work vector 371a7e14dcfSSatish Balay 372d8d19677SJose E. Roman Output Parameters: 373a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result 374a7e14dcfSSatish Balay . Db - row scaling component of the result 375a7e14dcfSSatish Balay - Dm - derivative with respect to scaling parameter 376a7e14dcfSSatish Balay 377a7e14dcfSSatish Balay Level: developer 378a7e14dcfSSatish Balay 379db781477SPatrick Sanan .seealso `MatDFischer()` 380a7e14dcfSSatish Balay @*/ 381*9371c9d4SSatish Balay PetscErrorCode MatDSFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, PetscReal mu, Vec T1, Vec T2, Vec Da, Vec Db, Vec Dm) { 382a7e14dcfSSatish Balay PetscInt i, nn; 38346bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 38446bdf8c8SLisandro Dalcin PetscScalar *da, *db, *dm; 385a7e14dcfSSatish Balay PetscReal ai, bi, ci, di, ei, fi; 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay PetscFunctionBegin; 388a7e14dcfSSatish Balay if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) { 3899566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Dm)); 3909566063dSJacob Faibussowitsch PetscCall(MatDFischer(jac, X, Con, XL, XU, T1, T2, Da, Db)); 3912d0e5244SBarry Smith } else { 3929566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nn)); 3939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 3949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Con, &f)); 3959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &l)); 3969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &u)); 3979566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 3989566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 3999566063dSJacob Faibussowitsch PetscCall(VecGetArray(Dm, &dm)); 400a7e14dcfSSatish Balay 401a7e14dcfSSatish Balay for (i = 0; i < nn; ++i) { 40246bdf8c8SLisandro Dalcin if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { 403a7e14dcfSSatish Balay da[i] = -mu; 40446bdf8c8SLisandro Dalcin db[i] = -1.0; 405a7e14dcfSSatish Balay dm[i] = -x[i]; 40646bdf8c8SLisandro Dalcin } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { 407658c1fc4SLisandro Dalcin bi = PetscRealPart(u[i]) - PetscRealPart(x[i]); 408658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 409a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 410a7e14dcfSSatish Balay 41146bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 412658c1fc4SLisandro Dalcin db[i] = -PetscRealPart(f[i]) / ai - 1.0; 413a7e14dcfSSatish Balay dm[i] = 2.0 * mu / ai; 41446bdf8c8SLisandro Dalcin } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { 415658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(l[i]); 416658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 417a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 418a7e14dcfSSatish Balay 41946bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 420658c1fc4SLisandro Dalcin db[i] = PetscRealPart(f[i]) / ai - 1.0; 421a7e14dcfSSatish Balay dm[i] = 2.0 * mu / ai; 422658c1fc4SLisandro Dalcin } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) { 42346bdf8c8SLisandro Dalcin da[i] = -1.0; 42446bdf8c8SLisandro Dalcin db[i] = 0.0; 42546bdf8c8SLisandro Dalcin dm[i] = 0.0; 4262d0e5244SBarry Smith } else { 427658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(u[i]); 428658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 429a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 430a7e14dcfSSatish Balay 43146bdf8c8SLisandro Dalcin ci = bi / ai + 1.0; 432658c1fc4SLisandro Dalcin di = PetscRealPart(f[i]) / ai + 1.0; 433a7e14dcfSSatish Balay fi = 2.0 * mu / ai; 434a7e14dcfSSatish Balay 435658c1fc4SLisandro Dalcin ei = SFischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]), mu); 436658c1fc4SLisandro Dalcin ai = fischsnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei, mu); 437a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 438a7e14dcfSSatish Balay 43946bdf8c8SLisandro Dalcin bi = ei / ai - 1.0; 440a7e14dcfSSatish Balay ei = 2.0 * mu / ei; 441658c1fc4SLisandro Dalcin ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0; 442a7e14dcfSSatish Balay 443a7e14dcfSSatish Balay da[i] = ai + bi * ci; 444a7e14dcfSSatish Balay db[i] = bi * di; 445a7e14dcfSSatish Balay dm[i] = ei + bi * fi; 446a7e14dcfSSatish Balay } 447a7e14dcfSSatish Balay } 448a7e14dcfSSatish Balay 4499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 4509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Con, &f)); 4519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &l)); 4529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &u)); 4539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 4549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 4559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Dm, &dm)); 456a7e14dcfSSatish Balay } 457a7e14dcfSSatish Balay PetscFunctionReturn(0); 458a7e14dcfSSatish Balay } 459a7e14dcfSSatish Balay 460*9371c9d4SSatish Balay static inline PetscReal ST_InternalPN(PetscScalar in, PetscReal lb, PetscReal ub) { 4618370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in) - ub) - PetscMax(0, -(PetscReal)PetscRealPart(in) - PetscAbsReal(lb)); 4628370d7cdSHansol Suh } 4638370d7cdSHansol Suh 464*9371c9d4SSatish Balay static inline PetscReal ST_InternalNN(PetscScalar in, PetscReal lb, PetscReal ub) { 4658370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in) + PetscAbsReal(ub)) - PetscMax(0, -(PetscReal)PetscRealPart(in) - PetscAbsReal(lb)); 4668370d7cdSHansol Suh } 4678370d7cdSHansol Suh 468*9371c9d4SSatish Balay static inline PetscReal ST_InternalPP(PetscScalar in, PetscReal lb, PetscReal ub) { 4698370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in) - ub) + PetscMin(0, (PetscReal)PetscRealPart(in) - lb); 4708370d7cdSHansol Suh } 4718370d7cdSHansol Suh 4728370d7cdSHansol Suh /*@ 4738370d7cdSHansol Suh TaoSoftThreshold - Calculates soft thresholding routine with input vector 4748370d7cdSHansol Suh and given lower and upper bound and returns it to output vector. 4758370d7cdSHansol Suh 4768370d7cdSHansol Suh Input Parameters: 4778370d7cdSHansol Suh + in - input vector to be thresholded 4788370d7cdSHansol Suh . lb - lower bound 479f0fc11ceSJed Brown - ub - upper bound 4808370d7cdSHansol Suh 48197bb3fdcSJose E. Roman Output Parameter: 4828370d7cdSHansol Suh . out - Soft thresholded output vector 4838370d7cdSHansol Suh 4848370d7cdSHansol Suh Notes: 4858370d7cdSHansol Suh Soft thresholding is defined as 4868370d7cdSHansol Suh \[ S(input,lb,ub) = 4878370d7cdSHansol Suh \begin{cases} 4888370d7cdSHansol Suh input - ub \text{input > ub} \\ 4898370d7cdSHansol Suh 0 \text{lb =< input <= ub} \\ 4908370d7cdSHansol Suh input + lb \text{input < lb} \\ 4918370d7cdSHansol Suh \] 4928370d7cdSHansol Suh 4938370d7cdSHansol Suh Level: developer 4948370d7cdSHansol Suh 4958370d7cdSHansol Suh @*/ 496*9371c9d4SSatish Balay PetscErrorCode TaoSoftThreshold(Vec in, PetscReal lb, PetscReal ub, Vec out) { 4978370d7cdSHansol Suh PetscInt i, nlocal, mlocal; 4988370d7cdSHansol Suh PetscScalar *inarray, *outarray; 4998370d7cdSHansol Suh 5008370d7cdSHansol Suh PetscFunctionBegin; 5019566063dSJacob Faibussowitsch PetscCall(VecGetArrayPair(in, out, &inarray, &outarray)); 5029566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(in, &nlocal)); 5039566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(in, &mlocal)); 5048370d7cdSHansol Suh 5053c859ba3SBarry Smith PetscCheck(nlocal == mlocal, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input and output vectors need to be of same size."); 5063c859ba3SBarry Smith PetscCheck(lb < ub, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Lower bound needs to be lower than upper bound."); 5078370d7cdSHansol Suh 5088370d7cdSHansol Suh if (ub >= 0 && lb < 0) { 5098370d7cdSHansol Suh for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalPN(inarray[i], lb, ub); 5108370d7cdSHansol Suh } else if (ub < 0 && lb < 0) { 5118370d7cdSHansol Suh for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalNN(inarray[i], lb, ub); 5128370d7cdSHansol Suh } else { 5138370d7cdSHansol Suh for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalPP(inarray[i], lb, ub); 5148370d7cdSHansol Suh } 5158370d7cdSHansol Suh 5169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayPair(in, out, &inarray, &outarray)); 5178370d7cdSHansol Suh PetscFunctionReturn(0); 5188370d7cdSHansol Suh } 519