1af0996ceSBarry Smith #include <petsc/private/petscimpl.h> 2ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 38370d7cdSHansol Suh #include <petscsys.h> 4a7e14dcfSSatish Balay 5d71ae5a4SJacob Faibussowitsch static inline PetscReal Fischer(PetscReal a, PetscReal b) 6d71ae5a4SJacob Faibussowitsch { 7a7e14dcfSSatish Balay /* Method suggested by Bob Vanderbei */ 8ad540459SPierre Jolivet if (a + b <= 0) return PetscSqrtReal(a * a + b * b) - (a + b); 946bdf8c8SLisandro Dalcin return -2.0 * a * b / (PetscSqrtReal(a * a + b * b) + (a + b)); 10a7e14dcfSSatish Balay } 11a7e14dcfSSatish Balay 12a7e14dcfSSatish Balay /*@ 13a7e14dcfSSatish Balay VecFischer - Evaluates the Fischer-Burmeister function for complementarity 14a7e14dcfSSatish Balay problems. 15a7e14dcfSSatish Balay 16*a1cb98faSBarry Smith Logically Collective on X 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay Input Parameters: 19a7e14dcfSSatish Balay + X - current point 20a7e14dcfSSatish Balay . F - function evaluated at x 21a7e14dcfSSatish Balay . L - lower bounds 22a7e14dcfSSatish Balay - U - upper bounds 23a7e14dcfSSatish Balay 24f899ff85SJose E. Roman Output Parameter: 25a7e14dcfSSatish Balay . FB - The Fischer-Burmeister function vector 26a7e14dcfSSatish Balay 27a7e14dcfSSatish Balay Notes: 28a7e14dcfSSatish Balay The Fischer-Burmeister function is defined as 29a7e14dcfSSatish Balay $ phi(a,b) := sqrt(a*a + b*b) - a - b 30a7e14dcfSSatish Balay and is used reformulate a complementarity problem as a semismooth 31a7e14dcfSSatish Balay system of equations. 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay The result of this function is done by cases: 34a7e14dcfSSatish Balay + l[i] == -infinity, u[i] == infinity -- fb[i] = -f[i] 35a7e14dcfSSatish Balay . l[i] == -infinity, u[i] finite -- fb[i] = phi(u[i]-x[i], -f[i]) 36a7e14dcfSSatish Balay . l[i] finite, u[i] == infinity -- fb[i] = phi(x[i]-l[i], f[i]) 37a7e14dcfSSatish Balay . l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u])) 38a7e14dcfSSatish Balay - otherwise l[i] == u[i] -- fb[i] = l[i] - x[i] 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay Level: developer 41a7e14dcfSSatish Balay 42*a1cb98faSBarry Smith .seealso: `Vec`, `VecSFischer()`, `MatDFischer()`, `MatDSFischer()` 43a7e14dcfSSatish Balay @*/ 44d71ae5a4SJacob Faibussowitsch PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB) 45d71ae5a4SJacob Faibussowitsch { 4646bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 4746bdf8c8SLisandro Dalcin PetscScalar *fb; 48a7e14dcfSSatish Balay PetscReal xval, fval, lval, uval; 49a7e14dcfSSatish Balay PetscInt low[5], high[5], n, i; 50a7e14dcfSSatish Balay 51a7e14dcfSSatish Balay PetscFunctionBegin; 52a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 53a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID, 2); 5476be6f4fSStefano Zampini if (L) PetscValidHeaderSpecific(L, VEC_CLASSID, 3); 5576be6f4fSStefano Zampini if (U) PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 56064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(FB, VEC_CLASSID, 5); 57a7e14dcfSSatish Balay 5876be6f4fSStefano Zampini if (!L && !U) { 5976be6f4fSStefano Zampini PetscCall(VecAXPBY(FB, -1.0, 0.0, F)); 6076be6f4fSStefano Zampini PetscFunctionReturn(0); 6176be6f4fSStefano Zampini } 6276be6f4fSStefano Zampini 639566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, low, high)); 649566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(F, low + 1, high + 1)); 659566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(L, low + 2, high + 2)); 669566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(U, low + 3, high + 3)); 679566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4)); 68a7e14dcfSSatish Balay 69ad540459SPierre Jolivet 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"); 70a7e14dcfSSatish Balay 719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(L, &l)); 749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 759566063dSJacob Faibussowitsch PetscCall(VecGetArray(FB, &fb)); 76a7e14dcfSSatish Balay 779566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay for (i = 0; i < n; ++i) { 8076be6f4fSStefano Zampini xval = PetscRealPart(x[i]); 8176be6f4fSStefano Zampini fval = PetscRealPart(f[i]); 8276be6f4fSStefano Zampini lval = PetscRealPart(l[i]); 8376be6f4fSStefano Zampini uval = PetscRealPart(u[i]); 84a7e14dcfSSatish Balay 8576be6f4fSStefano Zampini if (lval <= -PETSC_INFINITY && uval >= PETSC_INFINITY) { 86a7e14dcfSSatish Balay fb[i] = -fval; 87e270355aSBarry Smith } else if (lval <= -PETSC_INFINITY) { 88a7e14dcfSSatish Balay fb[i] = -Fischer(uval - xval, -fval); 89e270355aSBarry Smith } else if (uval >= PETSC_INFINITY) { 90a7e14dcfSSatish Balay fb[i] = Fischer(xval - lval, fval); 912d0e5244SBarry Smith } else if (lval == uval) { 92a7e14dcfSSatish Balay fb[i] = lval - xval; 932d0e5244SBarry Smith } else { 94a7e14dcfSSatish Balay fval = Fischer(uval - xval, -fval); 95a7e14dcfSSatish Balay fb[i] = Fischer(xval - lval, fval); 96a7e14dcfSSatish Balay } 97a7e14dcfSSatish Balay } 98a7e14dcfSSatish Balay 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 1019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(L, &l)); 1029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(FB, &fb)); 104a7e14dcfSSatish Balay PetscFunctionReturn(0); 105a7e14dcfSSatish Balay } 106a7e14dcfSSatish Balay 107d71ae5a4SJacob Faibussowitsch static inline PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c) 108d71ae5a4SJacob Faibussowitsch { 109a7e14dcfSSatish Balay /* Method suggested by Bob Vanderbei */ 110ad540459SPierre Jolivet if (a + b <= 0) return PetscSqrtReal(a * a + b * b + 2.0 * c * c) - (a + b); 1113f6ba705SLisandro Dalcin return 2.0 * (c * c - a * b) / (PetscSqrtReal(a * a + b * b + 2.0 * c * c) + (a + b)); 112a7e14dcfSSatish Balay } 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay /*@ 115a7e14dcfSSatish Balay VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for 116a7e14dcfSSatish Balay complementarity problems. 117a7e14dcfSSatish Balay 118*a1cb98faSBarry Smith Logically Collective on X 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay Input Parameters: 121a7e14dcfSSatish Balay + X - current point 122a7e14dcfSSatish Balay . F - function evaluated at x 123a7e14dcfSSatish Balay . L - lower bounds 124a7e14dcfSSatish Balay . U - upper bounds 125a7e14dcfSSatish Balay - mu - smoothing parameter 126a7e14dcfSSatish Balay 127f899ff85SJose E. Roman Output Parameter: 128a7e14dcfSSatish Balay . FB - The Smoothed Fischer-Burmeister function vector 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay Notes: 131a7e14dcfSSatish Balay The Smoothed Fischer-Burmeister function is defined as 132a7e14dcfSSatish Balay $ phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b 133a7e14dcfSSatish Balay and is used reformulate a complementarity problem as a semismooth 134a7e14dcfSSatish Balay system of equations. 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay The result of this function is done by cases: 137a7e14dcfSSatish Balay + l[i] == -infinity, u[i] == infinity -- fb[i] = -f[i] - 2*mu*x[i] 138a7e14dcfSSatish Balay . l[i] == -infinity, u[i] finite -- fb[i] = phi(u[i]-x[i], -f[i], mu) 139a7e14dcfSSatish Balay . l[i] finite, u[i] == infinity -- fb[i] = phi(x[i]-l[i], f[i], mu) 140a7e14dcfSSatish Balay . l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu) 141a7e14dcfSSatish Balay - otherwise l[i] == u[i] -- fb[i] = l[i] - x[i] 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay Level: developer 144a7e14dcfSSatish Balay 145*a1cb98faSBarry Smith .seealso: `Vec`, `VecFischer()`, `MatDFischer()`, `MatDSFischer()` 146a7e14dcfSSatish Balay @*/ 147d71ae5a4SJacob Faibussowitsch PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB) 148d71ae5a4SJacob Faibussowitsch { 14946bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 15046bdf8c8SLisandro Dalcin PetscScalar *fb; 151a7e14dcfSSatish Balay PetscReal xval, fval, lval, uval; 152a7e14dcfSSatish Balay PetscInt low[5], high[5], n, i; 153a7e14dcfSSatish Balay 154a7e14dcfSSatish Balay PetscFunctionBegin; 155a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 156a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID, 2); 157a7e14dcfSSatish Balay PetscValidHeaderSpecific(L, VEC_CLASSID, 3); 158a7e14dcfSSatish Balay PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 159a7e14dcfSSatish Balay PetscValidHeaderSpecific(FB, VEC_CLASSID, 6); 160a7e14dcfSSatish Balay 1619566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, low, high)); 1629566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(F, low + 1, high + 1)); 1639566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(L, low + 2, high + 2)); 1649566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(U, low + 3, high + 3)); 1659566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(FB, low + 4, high + 4)); 166a7e14dcfSSatish Balay 167ad540459SPierre Jolivet 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"); 168a7e14dcfSSatish Balay 1699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 1719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(L, &l)); 1729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1739566063dSJacob Faibussowitsch PetscCall(VecGetArray(FB, &fb)); 174a7e14dcfSSatish Balay 1759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 176a7e14dcfSSatish Balay 177a7e14dcfSSatish Balay for (i = 0; i < n; ++i) { 1789371c9d4SSatish Balay xval = PetscRealPart(*x++); 1799371c9d4SSatish Balay fval = PetscRealPart(*f++); 1809371c9d4SSatish Balay lval = PetscRealPart(*l++); 1819371c9d4SSatish Balay uval = PetscRealPart(*u++); 182a7e14dcfSSatish Balay 183e270355aSBarry Smith if ((lval <= -PETSC_INFINITY) && (uval >= PETSC_INFINITY)) { 184a7e14dcfSSatish Balay (*fb++) = -fval - mu * xval; 185e270355aSBarry Smith } else if (lval <= -PETSC_INFINITY) { 186a7e14dcfSSatish Balay (*fb++) = -SFischer(uval - xval, -fval, mu); 187e270355aSBarry Smith } else if (uval >= PETSC_INFINITY) { 188a7e14dcfSSatish Balay (*fb++) = SFischer(xval - lval, fval, mu); 1892d0e5244SBarry Smith } else if (lval == uval) { 190a7e14dcfSSatish Balay (*fb++) = lval - xval; 1912d0e5244SBarry Smith } else { 192a7e14dcfSSatish Balay fval = SFischer(uval - xval, -fval, mu); 193a7e14dcfSSatish Balay (*fb++) = SFischer(xval - lval, fval, mu); 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay } 1969371c9d4SSatish Balay x -= n; 1979371c9d4SSatish Balay f -= n; 1989371c9d4SSatish Balay l -= n; 1999371c9d4SSatish Balay u -= n; 2009371c9d4SSatish Balay fb -= n; 201a7e14dcfSSatish Balay 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(L, &l)); 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(FB, &fb)); 207a7e14dcfSSatish Balay PetscFunctionReturn(0); 208a7e14dcfSSatish Balay } 209a7e14dcfSSatish Balay 210d71ae5a4SJacob Faibussowitsch static inline PetscReal fischnorm(PetscReal a, PetscReal b) 211d71ae5a4SJacob Faibussowitsch { 212658c1fc4SLisandro Dalcin return PetscSqrtReal(a * a + b * b); 213a7e14dcfSSatish Balay } 214a7e14dcfSSatish Balay 215d71ae5a4SJacob Faibussowitsch static inline PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c) 216d71ae5a4SJacob Faibussowitsch { 217658c1fc4SLisandro Dalcin return PetscSqrtReal(a * a + b * b + 2.0 * c * c); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay /*@ 221235fd6e6SBarry Smith MatDFischer - 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 241*a1cb98faSBarry Smith .seealso: `Mat`, `VecFischer()`, `VecSFischer()`, `MatDSFischer()` 242a7e14dcfSSatish Balay @*/ 243d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db) 244d71ae5a4SJacob Faibussowitsch { 245a7e14dcfSSatish Balay PetscInt i, nn; 24646bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u, *t2; 24746bdf8c8SLisandro Dalcin PetscScalar *da, *db, *t1; 248a7e14dcfSSatish Balay PetscReal ai, bi, ci, di, ei; 249a7e14dcfSSatish Balay 250a7e14dcfSSatish Balay PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nn)); 2529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Con, &f)); 2549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &l)); 2559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &u)); 2569566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 2579566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArray(T1, &t1)); 2599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(T2, &t2)); 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay for (i = 0; i < nn; i++) { 26246bdf8c8SLisandro Dalcin da[i] = 0.0; 26346bdf8c8SLisandro Dalcin db[i] = 0.0; 26446bdf8c8SLisandro Dalcin t1[i] = 0.0; 265a7e14dcfSSatish Balay 26646bdf8c8SLisandro Dalcin if (PetscAbsScalar(f[i]) <= PETSC_MACHINE_EPSILON) { 26746bdf8c8SLisandro Dalcin if (PetscRealPart(l[i]) > PETSC_NINFINITY && PetscAbsScalar(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) { 26846bdf8c8SLisandro Dalcin t1[i] = 1.0; 26946bdf8c8SLisandro Dalcin da[i] = 1.0; 270a7e14dcfSSatish Balay } 271a7e14dcfSSatish Balay 27246bdf8c8SLisandro Dalcin if (PetscRealPart(u[i]) < PETSC_INFINITY && PetscAbsScalar(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) { 27346bdf8c8SLisandro Dalcin t1[i] = 1.0; 27446bdf8c8SLisandro Dalcin db[i] = 1.0; 275a7e14dcfSSatish Balay } 276a7e14dcfSSatish Balay } 277a7e14dcfSSatish Balay } 278a7e14dcfSSatish Balay 2799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(T1, &t1)); 2809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(T2, &t2)); 2819566063dSJacob Faibussowitsch PetscCall(MatMult(jac, T1, T2)); 2829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(T2, &t2)); 283a7e14dcfSSatish Balay 284a7e14dcfSSatish Balay for (i = 0; i < nn; i++) { 28546bdf8c8SLisandro Dalcin if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { 28646bdf8c8SLisandro Dalcin da[i] = 0.0; 28746bdf8c8SLisandro Dalcin db[i] = -1.0; 28846bdf8c8SLisandro Dalcin } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { 28946bdf8c8SLisandro Dalcin if (PetscRealPart(db[i]) >= 1) { 290658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 291a7e14dcfSSatish Balay 29246bdf8c8SLisandro Dalcin da[i] = -1.0 / ai - 1.0; 29346bdf8c8SLisandro Dalcin db[i] = -t2[i] / ai - 1.0; 2942d0e5244SBarry Smith } else { 295658c1fc4SLisandro Dalcin bi = PetscRealPart(u[i]) - PetscRealPart(x[i]); 296658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 297a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 298a7e14dcfSSatish Balay 29946bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 30046bdf8c8SLisandro Dalcin db[i] = -f[i] / ai - 1.0; 301a7e14dcfSSatish Balay } 30246bdf8c8SLisandro Dalcin } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { 30346bdf8c8SLisandro Dalcin if (PetscRealPart(da[i]) >= 1) { 304658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 305a7e14dcfSSatish Balay 30646bdf8c8SLisandro Dalcin da[i] = 1.0 / ai - 1.0; 30746bdf8c8SLisandro Dalcin db[i] = t2[i] / ai - 1.0; 3082d0e5244SBarry Smith } else { 309658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(l[i]); 310658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 311a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 312a7e14dcfSSatish Balay 31346bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 31446bdf8c8SLisandro Dalcin db[i] = f[i] / ai - 1.0; 315a7e14dcfSSatish Balay } 316658c1fc4SLisandro Dalcin } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) { 31746bdf8c8SLisandro Dalcin da[i] = -1.0; 31846bdf8c8SLisandro Dalcin db[i] = 0.0; 3192d0e5244SBarry Smith } else { 32046bdf8c8SLisandro Dalcin if (PetscRealPart(db[i]) >= 1) { 321658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 322a7e14dcfSSatish Balay 32346bdf8c8SLisandro Dalcin ci = 1.0 / ai + 1.0; 324658c1fc4SLisandro Dalcin di = PetscRealPart(t2[i]) / ai + 1.0; 3252d0e5244SBarry Smith } else { 326658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(u[i]); 327658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 328a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 329a7e14dcfSSatish Balay 33046bdf8c8SLisandro Dalcin ci = bi / ai + 1.0; 331658c1fc4SLisandro Dalcin di = PetscRealPart(f[i]) / ai + 1.0; 332a7e14dcfSSatish Balay } 333a7e14dcfSSatish Balay 33446bdf8c8SLisandro Dalcin if (PetscRealPart(da[i]) >= 1) { 335658c1fc4SLisandro Dalcin bi = ci + di * PetscRealPart(t2[i]); 336658c1fc4SLisandro Dalcin ai = fischnorm(1.0, bi); 337a7e14dcfSSatish Balay 33846bdf8c8SLisandro Dalcin bi = bi / ai - 1.0; 33946bdf8c8SLisandro Dalcin ai = 1.0 / ai - 1.0; 3402d0e5244SBarry Smith } else { 341658c1fc4SLisandro Dalcin ei = Fischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i])); 342658c1fc4SLisandro Dalcin ai = fischnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei); 343a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 344a7e14dcfSSatish Balay 34546bdf8c8SLisandro Dalcin bi = ei / ai - 1.0; 346658c1fc4SLisandro Dalcin ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0; 347a7e14dcfSSatish Balay } 348a7e14dcfSSatish Balay 349a7e14dcfSSatish Balay da[i] = ai + bi * ci; 350a7e14dcfSSatish Balay db[i] = bi * di; 351a7e14dcfSSatish Balay } 352a7e14dcfSSatish Balay } 353a7e14dcfSSatish Balay 3549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 3559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 3569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Con, &f)); 3589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &l)); 3599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &u)); 3609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(T2, &t2)); 361a7e14dcfSSatish Balay PetscFunctionReturn(0); 3628e3154b5SSatish Balay } 363a7e14dcfSSatish Balay 364a7e14dcfSSatish Balay /*@ 365235fd6e6SBarry Smith MatDSFischer - Calculates an element of the B-subdifferential of the 366a7e14dcfSSatish Balay smoothed Fischer-Burmeister function for complementarity problems. 367a7e14dcfSSatish Balay 368a7e14dcfSSatish Balay Collective on jac 369a7e14dcfSSatish Balay 370a7e14dcfSSatish Balay Input Parameters: 371a7e14dcfSSatish Balay + jac - the jacobian of f at X 372a7e14dcfSSatish Balay . X - current point 373a7e14dcfSSatish Balay . F - constraint function evaluated at X 374a7e14dcfSSatish Balay . XL - lower bounds 375a7e14dcfSSatish Balay . XU - upper bounds 376a7e14dcfSSatish Balay . mu - smoothing parameter 377a7e14dcfSSatish Balay . T1 - work vector 378a7e14dcfSSatish Balay - T2 - work vector 379a7e14dcfSSatish Balay 380d8d19677SJose E. Roman Output Parameters: 381a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result 382a7e14dcfSSatish Balay . Db - row scaling component of the result 383a7e14dcfSSatish Balay - Dm - derivative with respect to scaling parameter 384a7e14dcfSSatish Balay 385a7e14dcfSSatish Balay Level: developer 386a7e14dcfSSatish Balay 387*a1cb98faSBarry Smith .seealso: `Mat`, `VecFischer()`, `VecDFischer()`, `MatDFischer()` 388a7e14dcfSSatish Balay @*/ 389d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDSFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, PetscReal mu, Vec T1, Vec T2, Vec Da, Vec Db, Vec Dm) 390d71ae5a4SJacob Faibussowitsch { 391a7e14dcfSSatish Balay PetscInt i, nn; 39246bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 39346bdf8c8SLisandro Dalcin PetscScalar *da, *db, *dm; 394a7e14dcfSSatish Balay PetscReal ai, bi, ci, di, ei, fi; 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay PetscFunctionBegin; 397a7e14dcfSSatish Balay if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) { 3989566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Dm)); 3999566063dSJacob Faibussowitsch PetscCall(MatDFischer(jac, X, Con, XL, XU, T1, T2, Da, Db)); 4002d0e5244SBarry Smith } else { 4019566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nn)); 4029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 4039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Con, &f)); 4049566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &l)); 4059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &u)); 4069566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 4079566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 4089566063dSJacob Faibussowitsch PetscCall(VecGetArray(Dm, &dm)); 409a7e14dcfSSatish Balay 410a7e14dcfSSatish Balay for (i = 0; i < nn; ++i) { 41146bdf8c8SLisandro Dalcin if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { 412a7e14dcfSSatish Balay da[i] = -mu; 41346bdf8c8SLisandro Dalcin db[i] = -1.0; 414a7e14dcfSSatish Balay dm[i] = -x[i]; 41546bdf8c8SLisandro Dalcin } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { 416658c1fc4SLisandro Dalcin bi = PetscRealPart(u[i]) - PetscRealPart(x[i]); 417658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 418a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 419a7e14dcfSSatish Balay 42046bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 421658c1fc4SLisandro Dalcin db[i] = -PetscRealPart(f[i]) / ai - 1.0; 422a7e14dcfSSatish Balay dm[i] = 2.0 * mu / ai; 42346bdf8c8SLisandro Dalcin } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { 424658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(l[i]); 425658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 426a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 427a7e14dcfSSatish Balay 42846bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 429658c1fc4SLisandro Dalcin db[i] = PetscRealPart(f[i]) / ai - 1.0; 430a7e14dcfSSatish Balay dm[i] = 2.0 * mu / ai; 431658c1fc4SLisandro Dalcin } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) { 43246bdf8c8SLisandro Dalcin da[i] = -1.0; 43346bdf8c8SLisandro Dalcin db[i] = 0.0; 43446bdf8c8SLisandro Dalcin dm[i] = 0.0; 4352d0e5244SBarry Smith } else { 436658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(u[i]); 437658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 438a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 439a7e14dcfSSatish Balay 44046bdf8c8SLisandro Dalcin ci = bi / ai + 1.0; 441658c1fc4SLisandro Dalcin di = PetscRealPart(f[i]) / ai + 1.0; 442a7e14dcfSSatish Balay fi = 2.0 * mu / ai; 443a7e14dcfSSatish Balay 444658c1fc4SLisandro Dalcin ei = SFischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]), mu); 445658c1fc4SLisandro Dalcin ai = fischsnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei, mu); 446a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 447a7e14dcfSSatish Balay 44846bdf8c8SLisandro Dalcin bi = ei / ai - 1.0; 449a7e14dcfSSatish Balay ei = 2.0 * mu / ei; 450658c1fc4SLisandro Dalcin ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0; 451a7e14dcfSSatish Balay 452a7e14dcfSSatish Balay da[i] = ai + bi * ci; 453a7e14dcfSSatish Balay db[i] = bi * di; 454a7e14dcfSSatish Balay dm[i] = ei + bi * fi; 455a7e14dcfSSatish Balay } 456a7e14dcfSSatish Balay } 457a7e14dcfSSatish Balay 4589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 4599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Con, &f)); 4609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &l)); 4619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &u)); 4629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 4639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 4649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Dm, &dm)); 465a7e14dcfSSatish Balay } 466a7e14dcfSSatish Balay PetscFunctionReturn(0); 467a7e14dcfSSatish Balay } 468a7e14dcfSSatish Balay 469d71ae5a4SJacob Faibussowitsch static inline PetscReal ST_InternalPN(PetscScalar in, PetscReal lb, PetscReal ub) 470d71ae5a4SJacob Faibussowitsch { 4718370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in) - ub) - PetscMax(0, -(PetscReal)PetscRealPart(in) - PetscAbsReal(lb)); 4728370d7cdSHansol Suh } 4738370d7cdSHansol Suh 474d71ae5a4SJacob Faibussowitsch static inline PetscReal ST_InternalNN(PetscScalar in, PetscReal lb, PetscReal ub) 475d71ae5a4SJacob Faibussowitsch { 4768370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in) + PetscAbsReal(ub)) - PetscMax(0, -(PetscReal)PetscRealPart(in) - PetscAbsReal(lb)); 4778370d7cdSHansol Suh } 4788370d7cdSHansol Suh 479d71ae5a4SJacob Faibussowitsch static inline PetscReal ST_InternalPP(PetscScalar in, PetscReal lb, PetscReal ub) 480d71ae5a4SJacob Faibussowitsch { 4818370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in) - ub) + PetscMin(0, (PetscReal)PetscRealPart(in) - lb); 4828370d7cdSHansol Suh } 4838370d7cdSHansol Suh 4848370d7cdSHansol Suh /*@ 4858370d7cdSHansol Suh TaoSoftThreshold - Calculates soft thresholding routine with input vector 4868370d7cdSHansol Suh and given lower and upper bound and returns it to output vector. 4878370d7cdSHansol Suh 4888370d7cdSHansol Suh Input Parameters: 4898370d7cdSHansol Suh + in - input vector to be thresholded 4908370d7cdSHansol Suh . lb - lower bound 491f0fc11ceSJed Brown - ub - upper bound 4928370d7cdSHansol Suh 49397bb3fdcSJose E. Roman Output Parameter: 4948370d7cdSHansol Suh . out - Soft thresholded output vector 4958370d7cdSHansol Suh 4968370d7cdSHansol Suh Notes: 4978370d7cdSHansol Suh Soft thresholding is defined as 4988370d7cdSHansol Suh \[ S(input,lb,ub) = 4998370d7cdSHansol Suh \begin{cases} 5008370d7cdSHansol Suh input - ub \text{input > ub} \\ 5018370d7cdSHansol Suh 0 \text{lb =< input <= ub} \\ 5028370d7cdSHansol Suh input + lb \text{input < lb} \\ 5038370d7cdSHansol Suh \] 5048370d7cdSHansol Suh 5058370d7cdSHansol Suh Level: developer 5068370d7cdSHansol Suh 507*a1cb98faSBarry Smith .seealso: `Tao`, `Vec` 5088370d7cdSHansol Suh @*/ 509d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSoftThreshold(Vec in, PetscReal lb, PetscReal ub, Vec out) 510d71ae5a4SJacob Faibussowitsch { 5118370d7cdSHansol Suh PetscInt i, nlocal, mlocal; 5128370d7cdSHansol Suh PetscScalar *inarray, *outarray; 5138370d7cdSHansol Suh 5148370d7cdSHansol Suh PetscFunctionBegin; 5159566063dSJacob Faibussowitsch PetscCall(VecGetArrayPair(in, out, &inarray, &outarray)); 5169566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(in, &nlocal)); 5179566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(in, &mlocal)); 5188370d7cdSHansol Suh 5193c859ba3SBarry Smith PetscCheck(nlocal == mlocal, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input and output vectors need to be of same size."); 5203c859ba3SBarry Smith PetscCheck(lb < ub, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Lower bound needs to be lower than upper bound."); 5218370d7cdSHansol Suh 5228370d7cdSHansol Suh if (ub >= 0 && lb < 0) { 5238370d7cdSHansol Suh for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalPN(inarray[i], lb, ub); 5248370d7cdSHansol Suh } else if (ub < 0 && lb < 0) { 5258370d7cdSHansol Suh for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalNN(inarray[i], lb, ub); 5268370d7cdSHansol Suh } else { 5278370d7cdSHansol Suh for (i = 0; i < nlocal; i++) outarray[i] = ST_InternalPP(inarray[i], lb, ub); 5288370d7cdSHansol Suh } 5298370d7cdSHansol Suh 5309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayPair(in, out, &inarray, &outarray)); 5318370d7cdSHansol Suh PetscFunctionReturn(0); 5328370d7cdSHansol Suh } 533