1af0996ceSBarry Smith #include <petsc/private/petscimpl.h> 2ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 38370d7cdSHansol Suh #include <petscsys.h> 4a7e14dcfSSatish Balay 59fbee547SJacob Faibussowitsch static inline PetscReal Fischer(PetscReal a, PetscReal b) 6a7e14dcfSSatish Balay { 7a7e14dcfSSatish Balay /* Method suggested by Bob Vanderbei */ 8a7e14dcfSSatish Balay if (a + b <= 0) { 946bdf8c8SLisandro Dalcin return PetscSqrtReal(a*a + b*b) - (a + b); 10a7e14dcfSSatish Balay } 1146bdf8c8SLisandro Dalcin return -2.0*a*b / (PetscSqrtReal(a*a + b*b) + (a + b)); 12a7e14dcfSSatish Balay } 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay /*@ 15a7e14dcfSSatish Balay VecFischer - Evaluates the Fischer-Burmeister function for complementarity 16a7e14dcfSSatish Balay problems. 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay Logically Collective on vectors 19a7e14dcfSSatish Balay 20a7e14dcfSSatish Balay Input Parameters: 21a7e14dcfSSatish Balay + X - current point 22a7e14dcfSSatish Balay . F - function evaluated at x 23a7e14dcfSSatish Balay . L - lower bounds 24a7e14dcfSSatish Balay - U - upper bounds 25a7e14dcfSSatish Balay 26f899ff85SJose E. Roman Output Parameter: 27a7e14dcfSSatish Balay . FB - The Fischer-Burmeister function vector 28a7e14dcfSSatish Balay 29a7e14dcfSSatish Balay Notes: 30a7e14dcfSSatish Balay The Fischer-Burmeister function is defined as 31a7e14dcfSSatish Balay $ phi(a,b) := sqrt(a*a + b*b) - a - b 32a7e14dcfSSatish Balay and is used reformulate a complementarity problem as a semismooth 33a7e14dcfSSatish Balay system of equations. 34a7e14dcfSSatish Balay 35a7e14dcfSSatish Balay The result of this function is done by cases: 36a7e14dcfSSatish Balay + l[i] == -infinity, u[i] == infinity -- fb[i] = -f[i] 37a7e14dcfSSatish Balay . l[i] == -infinity, u[i] finite -- fb[i] = phi(u[i]-x[i], -f[i]) 38a7e14dcfSSatish Balay . l[i] finite, u[i] == infinity -- fb[i] = phi(x[i]-l[i], f[i]) 39a7e14dcfSSatish Balay . l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u])) 40a7e14dcfSSatish Balay - otherwise l[i] == u[i] -- fb[i] = l[i] - x[i] 41a7e14dcfSSatish Balay 42a7e14dcfSSatish Balay Level: developer 43a7e14dcfSSatish Balay 44a7e14dcfSSatish Balay @*/ 45a7e14dcfSSatish Balay PetscErrorCode VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FB) 46a7e14dcfSSatish Balay { 4746bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 4846bdf8c8SLisandro Dalcin PetscScalar *fb; 49a7e14dcfSSatish Balay PetscReal xval, fval, lval, uval; 50a7e14dcfSSatish Balay PetscErrorCode ierr; 51a7e14dcfSSatish Balay PetscInt low[5], high[5], n, i; 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay PetscFunctionBegin; 54a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID,1); 55a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID,2); 56a7e14dcfSSatish Balay PetscValidHeaderSpecific(L, VEC_CLASSID,3); 57a7e14dcfSSatish Balay PetscValidHeaderSpecific(U, VEC_CLASSID,4); 58064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(FB, VEC_CLASSID,5); 59a7e14dcfSSatish Balay 60a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr); 61a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr); 62a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr); 65a7e14dcfSSatish Balay 66a7e14dcfSSatish Balay for (i = 1; i < 4; ++i) { 67*3c859ba3SBarry Smith PetscCheck(low[0] == low[i] && high[0] == high[i],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vectors must be identically loaded over processors"); 68a7e14dcfSSatish Balay } 69a7e14dcfSSatish Balay 705e081366SBarry Smith ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 715e081366SBarry Smith ierr = VecGetArrayRead(F, &f);CHKERRQ(ierr); 725e081366SBarry Smith ierr = VecGetArrayRead(L, &l);CHKERRQ(ierr); 735e081366SBarry Smith ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = VecGetArray(FB, &fb);CHKERRQ(ierr); 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr); 77a7e14dcfSSatish Balay 78a7e14dcfSSatish Balay for (i = 0; i < n; ++i) { 79658c1fc4SLisandro Dalcin xval = PetscRealPart(x[i]); fval = PetscRealPart(f[i]); 80658c1fc4SLisandro Dalcin lval = PetscRealPart(l[i]); uval = PetscRealPart(u[i]); 81a7e14dcfSSatish Balay 82e270355aSBarry Smith 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 965e081366SBarry Smith ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 975e081366SBarry Smith ierr = VecRestoreArrayRead(F, &f);CHKERRQ(ierr); 985e081366SBarry Smith ierr = VecRestoreArrayRead(L, &l);CHKERRQ(ierr); 995e081366SBarry Smith ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 100a7e14dcfSSatish Balay ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr); 101a7e14dcfSSatish Balay PetscFunctionReturn(0); 102a7e14dcfSSatish Balay } 103a7e14dcfSSatish Balay 1049fbee547SJacob Faibussowitsch static inline PetscReal SFischer(PetscReal a, PetscReal b, PetscReal c) 105a7e14dcfSSatish Balay { 106a7e14dcfSSatish Balay /* Method suggested by Bob Vanderbei */ 107a7e14dcfSSatish Balay if (a + b <= 0) { 1083f6ba705SLisandro Dalcin return PetscSqrtReal(a*a + b*b + 2.0*c*c) - (a + b); 109a7e14dcfSSatish Balay } 1103f6ba705SLisandro Dalcin return 2.0*(c*c - a*b) / (PetscSqrtReal(a*a + b*b + 2.0*c*c) + (a + b)); 111a7e14dcfSSatish Balay } 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay /*@ 114a7e14dcfSSatish Balay VecSFischer - Evaluates the Smoothed Fischer-Burmeister function for 115a7e14dcfSSatish Balay complementarity problems. 116a7e14dcfSSatish Balay 117a7e14dcfSSatish Balay Logically Collective on vectors 118a7e14dcfSSatish Balay 119a7e14dcfSSatish Balay Input Parameters: 120a7e14dcfSSatish Balay + X - current point 121a7e14dcfSSatish Balay . F - function evaluated at x 122a7e14dcfSSatish Balay . L - lower bounds 123a7e14dcfSSatish Balay . U - upper bounds 124a7e14dcfSSatish Balay - mu - smoothing parameter 125a7e14dcfSSatish Balay 126f899ff85SJose E. Roman Output Parameter: 127a7e14dcfSSatish Balay . FB - The Smoothed Fischer-Burmeister function vector 128a7e14dcfSSatish Balay 129a7e14dcfSSatish Balay Notes: 130a7e14dcfSSatish Balay The Smoothed Fischer-Burmeister function is defined as 131a7e14dcfSSatish Balay $ phi(a,b) := sqrt(a*a + b*b + 2*mu*mu) - a - b 132a7e14dcfSSatish Balay and is used reformulate a complementarity problem as a semismooth 133a7e14dcfSSatish Balay system of equations. 134a7e14dcfSSatish Balay 135a7e14dcfSSatish Balay The result of this function is done by cases: 136a7e14dcfSSatish Balay + l[i] == -infinity, u[i] == infinity -- fb[i] = -f[i] - 2*mu*x[i] 137a7e14dcfSSatish Balay . l[i] == -infinity, u[i] finite -- fb[i] = phi(u[i]-x[i], -f[i], mu) 138a7e14dcfSSatish Balay . l[i] finite, u[i] == infinity -- fb[i] = phi(x[i]-l[i], f[i], mu) 139a7e14dcfSSatish Balay . l[i] finite < u[i] finite -- fb[i] = phi(x[i]-l[i], phi(u[i]-x[i], -f[u], mu), mu) 140a7e14dcfSSatish Balay - otherwise l[i] == u[i] -- fb[i] = l[i] - x[i] 141a7e14dcfSSatish Balay 142a7e14dcfSSatish Balay Level: developer 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay .seealso VecFischer() 145a7e14dcfSSatish Balay @*/ 146a7e14dcfSSatish Balay PetscErrorCode VecSFischer(Vec X, Vec F, Vec L, Vec U, PetscReal mu, Vec FB) 147a7e14dcfSSatish Balay { 14846bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 14946bdf8c8SLisandro Dalcin PetscScalar *fb; 150a7e14dcfSSatish Balay PetscReal xval, fval, lval, uval; 151a7e14dcfSSatish Balay PetscErrorCode ierr; 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 161a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(X, low, high);CHKERRQ(ierr); 162a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(F, low + 1, high + 1);CHKERRQ(ierr); 163a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(L, low + 2, high + 2);CHKERRQ(ierr); 164a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(U, low + 3, high + 3);CHKERRQ(ierr); 165a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(FB, low + 4, high + 4);CHKERRQ(ierr); 166a7e14dcfSSatish Balay 167a7e14dcfSSatish Balay for (i = 1; i < 4; ++i) { 168*3c859ba3SBarry Smith PetscCheck(low[0] == low[i] && high[0] == high[i],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vectors must be identically loaded over processors"); 169a7e14dcfSSatish Balay } 170a7e14dcfSSatish Balay 1715e081366SBarry Smith ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 1725e081366SBarry Smith ierr = VecGetArrayRead(F, &f);CHKERRQ(ierr); 1735e081366SBarry Smith ierr = VecGetArrayRead(L, &l);CHKERRQ(ierr); 1745e081366SBarry Smith ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 175a7e14dcfSSatish Balay ierr = VecGetArray(FB, &fb);CHKERRQ(ierr); 176a7e14dcfSSatish Balay 177a7e14dcfSSatish Balay ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr); 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay for (i = 0; i < n; ++i) { 180658c1fc4SLisandro Dalcin xval = PetscRealPart(*x++); fval = PetscRealPart(*f++); 181658c1fc4SLisandro Dalcin lval = PetscRealPart(*l++); 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 } 196a7e14dcfSSatish Balay x -= n; f -= n; l -=n; u -= n; fb -= n; 197a7e14dcfSSatish Balay 1985e081366SBarry Smith ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 1995e081366SBarry Smith ierr = VecRestoreArrayRead(F, &f);CHKERRQ(ierr); 2005e081366SBarry Smith ierr = VecRestoreArrayRead(L, &l);CHKERRQ(ierr); 2015e081366SBarry Smith ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 202a7e14dcfSSatish Balay ierr = VecRestoreArray(FB, &fb);CHKERRQ(ierr); 203a7e14dcfSSatish Balay PetscFunctionReturn(0); 204a7e14dcfSSatish Balay } 205a7e14dcfSSatish Balay 2069fbee547SJacob Faibussowitsch static inline PetscReal fischnorm(PetscReal a, PetscReal b) 207a7e14dcfSSatish Balay { 208658c1fc4SLisandro Dalcin return PetscSqrtReal(a*a + b*b); 209a7e14dcfSSatish Balay } 210a7e14dcfSSatish Balay 2119fbee547SJacob Faibussowitsch static inline PetscReal fischsnorm(PetscReal a, PetscReal b, PetscReal c) 212a7e14dcfSSatish Balay { 213658c1fc4SLisandro Dalcin return PetscSqrtReal(a*a + b*b + 2.0*c*c); 214a7e14dcfSSatish Balay } 215a7e14dcfSSatish Balay 216a7e14dcfSSatish Balay /*@ 217235fd6e6SBarry Smith MatDFischer - Calculates an element of the B-subdifferential of the 218a7e14dcfSSatish Balay Fischer-Burmeister function for complementarity problems. 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay Collective on jac 221a7e14dcfSSatish Balay 222a7e14dcfSSatish Balay Input Parameters: 223a7e14dcfSSatish Balay + jac - the jacobian of f at X 224a7e14dcfSSatish Balay . X - current point 225a7e14dcfSSatish Balay . Con - constraints function evaluated at X 226a7e14dcfSSatish Balay . XL - lower bounds 227a7e14dcfSSatish Balay . XU - upper bounds 228a7e14dcfSSatish Balay . t1 - work vector 229a7e14dcfSSatish Balay - t2 - work vector 230a7e14dcfSSatish Balay 231a7e14dcfSSatish Balay Output Parameters: 232a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result 233a7e14dcfSSatish Balay - Db - row scaling component of the result 234a7e14dcfSSatish Balay 235a7e14dcfSSatish Balay Level: developer 236a7e14dcfSSatish Balay 237a7e14dcfSSatish Balay .seealso: VecFischer() 238a7e14dcfSSatish Balay @*/ 239235fd6e6SBarry Smith PetscErrorCode MatDFischer(Mat jac, Vec X, Vec Con, Vec XL, Vec XU, Vec T1, Vec T2, Vec Da, Vec Db) 240a7e14dcfSSatish Balay { 241a7e14dcfSSatish Balay PetscErrorCode ierr; 242a7e14dcfSSatish Balay PetscInt i,nn; 24346bdf8c8SLisandro Dalcin const PetscScalar *x,*f,*l,*u,*t2; 24446bdf8c8SLisandro Dalcin PetscScalar *da,*db,*t1; 245a7e14dcfSSatish Balay PetscReal ai,bi,ci,di,ei; 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay PetscFunctionBegin; 248a7e14dcfSSatish Balay ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr); 2495e081366SBarry Smith ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 2505e081366SBarry Smith ierr = VecGetArrayRead(Con,&f);CHKERRQ(ierr); 2515e081366SBarry Smith ierr = VecGetArrayRead(XL,&l);CHKERRQ(ierr); 2525e081366SBarry Smith ierr = VecGetArrayRead(XU,&u);CHKERRQ(ierr); 253a7e14dcfSSatish Balay ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 254a7e14dcfSSatish Balay ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 255a7e14dcfSSatish Balay ierr = VecGetArray(T1,&t1);CHKERRQ(ierr); 2565e081366SBarry Smith ierr = VecGetArrayRead(T2,&t2);CHKERRQ(ierr); 257a7e14dcfSSatish Balay 258a7e14dcfSSatish Balay for (i = 0; i < nn; i++) { 25946bdf8c8SLisandro Dalcin da[i] = 0.0; 26046bdf8c8SLisandro Dalcin db[i] = 0.0; 26146bdf8c8SLisandro Dalcin t1[i] = 0.0; 262a7e14dcfSSatish Balay 26346bdf8c8SLisandro Dalcin if (PetscAbsScalar(f[i]) <= PETSC_MACHINE_EPSILON) { 26446bdf8c8SLisandro Dalcin if (PetscRealPart(l[i]) > PETSC_NINFINITY && PetscAbsScalar(x[i] - l[i]) <= PETSC_MACHINE_EPSILON) { 26546bdf8c8SLisandro Dalcin t1[i] = 1.0; 26646bdf8c8SLisandro Dalcin da[i] = 1.0; 267a7e14dcfSSatish Balay } 268a7e14dcfSSatish Balay 26946bdf8c8SLisandro Dalcin if (PetscRealPart(u[i]) < PETSC_INFINITY && PetscAbsScalar(u[i] - x[i]) <= PETSC_MACHINE_EPSILON) { 27046bdf8c8SLisandro Dalcin t1[i] = 1.0; 27146bdf8c8SLisandro Dalcin db[i] = 1.0; 272a7e14dcfSSatish Balay } 273a7e14dcfSSatish Balay } 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay ierr = VecRestoreArray(T1,&t1);CHKERRQ(ierr); 2775e081366SBarry Smith ierr = VecRestoreArrayRead(T2,&t2);CHKERRQ(ierr); 278a7e14dcfSSatish Balay ierr = MatMult(jac,T1,T2);CHKERRQ(ierr); 2795e081366SBarry Smith ierr = VecGetArrayRead(T2,&t2);CHKERRQ(ierr); 280a7e14dcfSSatish Balay 281a7e14dcfSSatish Balay for (i = 0; i < nn; i++) { 28246bdf8c8SLisandro Dalcin if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { 28346bdf8c8SLisandro Dalcin da[i] = 0.0; 28446bdf8c8SLisandro Dalcin db[i] = -1.0; 28546bdf8c8SLisandro Dalcin } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { 28646bdf8c8SLisandro Dalcin if (PetscRealPart(db[i]) >= 1) { 287658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 288a7e14dcfSSatish Balay 28946bdf8c8SLisandro Dalcin da[i] = -1.0 / ai - 1.0; 29046bdf8c8SLisandro Dalcin db[i] = -t2[i] / ai - 1.0; 2912d0e5244SBarry Smith } else { 292658c1fc4SLisandro Dalcin bi = PetscRealPart(u[i]) - PetscRealPart(x[i]); 293658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 294a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 295a7e14dcfSSatish Balay 29646bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 29746bdf8c8SLisandro Dalcin db[i] = -f[i] / ai - 1.0; 298a7e14dcfSSatish Balay } 29946bdf8c8SLisandro Dalcin } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { 30046bdf8c8SLisandro Dalcin if (PetscRealPart(da[i]) >= 1) { 301658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 302a7e14dcfSSatish Balay 30346bdf8c8SLisandro Dalcin da[i] = 1.0 / ai - 1.0; 30446bdf8c8SLisandro Dalcin db[i] = t2[i] / ai - 1.0; 3052d0e5244SBarry Smith } else { 306658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(l[i]); 307658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 308a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 309a7e14dcfSSatish Balay 31046bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 31146bdf8c8SLisandro Dalcin db[i] = f[i] / ai - 1.0; 312a7e14dcfSSatish Balay } 313658c1fc4SLisandro Dalcin } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) { 31446bdf8c8SLisandro Dalcin da[i] = -1.0; 31546bdf8c8SLisandro Dalcin db[i] = 0.0; 3162d0e5244SBarry Smith } else { 31746bdf8c8SLisandro Dalcin if (PetscRealPart(db[i]) >= 1) { 318658c1fc4SLisandro Dalcin ai = fischnorm(1.0, PetscRealPart(t2[i])); 319a7e14dcfSSatish Balay 32046bdf8c8SLisandro Dalcin ci = 1.0 / ai + 1.0; 321658c1fc4SLisandro Dalcin di = PetscRealPart(t2[i]) / ai + 1.0; 3222d0e5244SBarry Smith } else { 323658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(u[i]); 324658c1fc4SLisandro Dalcin ai = fischnorm(bi, PetscRealPart(f[i])); 325a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 326a7e14dcfSSatish Balay 32746bdf8c8SLisandro Dalcin ci = bi / ai + 1.0; 328658c1fc4SLisandro Dalcin di = PetscRealPart(f[i]) / ai + 1.0; 329a7e14dcfSSatish Balay } 330a7e14dcfSSatish Balay 33146bdf8c8SLisandro Dalcin if (PetscRealPart(da[i]) >= 1) { 332658c1fc4SLisandro Dalcin bi = ci + di*PetscRealPart(t2[i]); 333658c1fc4SLisandro Dalcin ai = fischnorm(1.0, bi); 334a7e14dcfSSatish Balay 33546bdf8c8SLisandro Dalcin bi = bi / ai - 1.0; 33646bdf8c8SLisandro Dalcin ai = 1.0 / ai - 1.0; 3372d0e5244SBarry Smith } else { 338658c1fc4SLisandro Dalcin ei = Fischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i])); 339658c1fc4SLisandro Dalcin ai = fischnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei); 340a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 341a7e14dcfSSatish Balay 34246bdf8c8SLisandro Dalcin bi = ei / ai - 1.0; 343658c1fc4SLisandro Dalcin ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0; 344a7e14dcfSSatish Balay } 345a7e14dcfSSatish Balay 346a7e14dcfSSatish Balay da[i] = ai + bi*ci; 347a7e14dcfSSatish Balay db[i] = bi*di; 348a7e14dcfSSatish Balay } 349a7e14dcfSSatish Balay } 350a7e14dcfSSatish Balay 351a7e14dcfSSatish Balay ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 352a7e14dcfSSatish Balay ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 3535e081366SBarry Smith ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 3545e081366SBarry Smith ierr = VecRestoreArrayRead(Con,&f);CHKERRQ(ierr); 3555e081366SBarry Smith ierr = VecRestoreArrayRead(XL,&l);CHKERRQ(ierr); 3565e081366SBarry Smith ierr = VecRestoreArrayRead(XU,&u);CHKERRQ(ierr); 3575e081366SBarry Smith ierr = VecRestoreArrayRead(T2,&t2);CHKERRQ(ierr); 358a7e14dcfSSatish Balay PetscFunctionReturn(0); 3598e3154b5SSatish Balay } 360a7e14dcfSSatish Balay 361a7e14dcfSSatish Balay /*@ 362235fd6e6SBarry Smith MatDSFischer - Calculates an element of the B-subdifferential of the 363a7e14dcfSSatish Balay smoothed Fischer-Burmeister function for complementarity problems. 364a7e14dcfSSatish Balay 365a7e14dcfSSatish Balay Collective on jac 366a7e14dcfSSatish Balay 367a7e14dcfSSatish Balay Input Parameters: 368a7e14dcfSSatish Balay + jac - the jacobian of f at X 369a7e14dcfSSatish Balay . X - current point 370a7e14dcfSSatish Balay . F - constraint function evaluated at X 371a7e14dcfSSatish Balay . XL - lower bounds 372a7e14dcfSSatish Balay . XU - upper bounds 373a7e14dcfSSatish Balay . mu - smoothing parameter 374a7e14dcfSSatish Balay . T1 - work vector 375a7e14dcfSSatish Balay - T2 - work vector 376a7e14dcfSSatish Balay 377d8d19677SJose E. Roman Output Parameters: 378a7e14dcfSSatish Balay + Da - diagonal perturbation component of the result 379a7e14dcfSSatish Balay . Db - row scaling component of the result 380a7e14dcfSSatish Balay - Dm - derivative with respect to scaling parameter 381a7e14dcfSSatish Balay 382a7e14dcfSSatish Balay Level: developer 383a7e14dcfSSatish Balay 384235fd6e6SBarry Smith .seealso MatDFischer() 385a7e14dcfSSatish Balay @*/ 386235fd6e6SBarry 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) 387a7e14dcfSSatish Balay { 388a7e14dcfSSatish Balay PetscErrorCode ierr; 389a7e14dcfSSatish Balay PetscInt i,nn; 39046bdf8c8SLisandro Dalcin const PetscScalar *x, *f, *l, *u; 39146bdf8c8SLisandro Dalcin PetscScalar *da, *db, *dm; 392a7e14dcfSSatish Balay PetscReal ai, bi, ci, di, ei, fi; 393a7e14dcfSSatish Balay 394a7e14dcfSSatish Balay PetscFunctionBegin; 395a7e14dcfSSatish Balay if (PetscAbsReal(mu) <= PETSC_MACHINE_EPSILON) { 396a7e14dcfSSatish Balay ierr = VecZeroEntries(Dm);CHKERRQ(ierr); 397235fd6e6SBarry Smith ierr = MatDFischer(jac, X, Con, XL, XU, T1, T2, Da, Db);CHKERRQ(ierr); 3982d0e5244SBarry Smith } else { 399a7e14dcfSSatish Balay ierr = VecGetLocalSize(X,&nn);CHKERRQ(ierr); 4005e081366SBarry Smith ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 4015e081366SBarry Smith ierr = VecGetArrayRead(Con,&f);CHKERRQ(ierr); 4025e081366SBarry Smith ierr = VecGetArrayRead(XL,&l);CHKERRQ(ierr); 4035e081366SBarry Smith ierr = VecGetArrayRead(XU,&u);CHKERRQ(ierr); 404a7e14dcfSSatish Balay ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 405a7e14dcfSSatish Balay ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 406a7e14dcfSSatish Balay ierr = VecGetArray(Dm,&dm);CHKERRQ(ierr); 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay for (i = 0; i < nn; ++i) { 40946bdf8c8SLisandro Dalcin if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { 410a7e14dcfSSatish Balay da[i] = -mu; 41146bdf8c8SLisandro Dalcin db[i] = -1.0; 412a7e14dcfSSatish Balay dm[i] = -x[i]; 41346bdf8c8SLisandro Dalcin } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { 414658c1fc4SLisandro Dalcin bi = PetscRealPart(u[i]) - PetscRealPart(x[i]); 415658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 416a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 417a7e14dcfSSatish Balay 41846bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 419658c1fc4SLisandro Dalcin db[i] = -PetscRealPart(f[i]) / ai - 1.0; 420a7e14dcfSSatish Balay dm[i] = 2.0 * mu / ai; 42146bdf8c8SLisandro Dalcin } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { 422658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(l[i]); 423658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 424a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 425a7e14dcfSSatish Balay 42646bdf8c8SLisandro Dalcin da[i] = bi / ai - 1.0; 427658c1fc4SLisandro Dalcin db[i] = PetscRealPart(f[i]) / ai - 1.0; 428a7e14dcfSSatish Balay dm[i] = 2.0 * mu / ai; 429658c1fc4SLisandro Dalcin } else if (PetscRealPart(l[i]) == PetscRealPart(u[i])) { 43046bdf8c8SLisandro Dalcin da[i] = -1.0; 43146bdf8c8SLisandro Dalcin db[i] = 0.0; 43246bdf8c8SLisandro Dalcin dm[i] = 0.0; 4332d0e5244SBarry Smith } else { 434658c1fc4SLisandro Dalcin bi = PetscRealPart(x[i]) - PetscRealPart(u[i]); 435658c1fc4SLisandro Dalcin ai = fischsnorm(bi, PetscRealPart(f[i]), mu); 436a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 437a7e14dcfSSatish Balay 43846bdf8c8SLisandro Dalcin ci = bi / ai + 1.0; 439658c1fc4SLisandro Dalcin di = PetscRealPart(f[i]) / ai + 1.0; 440a7e14dcfSSatish Balay fi = 2.0 * mu / ai; 441a7e14dcfSSatish Balay 442658c1fc4SLisandro Dalcin ei = SFischer(PetscRealPart(u[i]) - PetscRealPart(x[i]), -PetscRealPart(f[i]), mu); 443658c1fc4SLisandro Dalcin ai = fischsnorm(PetscRealPart(x[i]) - PetscRealPart(l[i]), ei, mu); 444a7e14dcfSSatish Balay ai = PetscMax(PETSC_MACHINE_EPSILON, ai); 445a7e14dcfSSatish Balay 44646bdf8c8SLisandro Dalcin bi = ei / ai - 1.0; 447a7e14dcfSSatish Balay ei = 2.0 * mu / ei; 448658c1fc4SLisandro Dalcin ai = (PetscRealPart(x[i]) - PetscRealPart(l[i])) / ai - 1.0; 449a7e14dcfSSatish Balay 450a7e14dcfSSatish Balay da[i] = ai + bi*ci; 451a7e14dcfSSatish Balay db[i] = bi*di; 452a7e14dcfSSatish Balay dm[i] = ei + bi*fi; 453a7e14dcfSSatish Balay } 454a7e14dcfSSatish Balay } 455a7e14dcfSSatish Balay 4565e081366SBarry Smith ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 4575e081366SBarry Smith ierr = VecRestoreArrayRead(Con,&f);CHKERRQ(ierr); 4585e081366SBarry Smith ierr = VecRestoreArrayRead(XL,&l);CHKERRQ(ierr); 4595e081366SBarry Smith ierr = VecRestoreArrayRead(XU,&u);CHKERRQ(ierr); 460a7e14dcfSSatish Balay ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 461a7e14dcfSSatish Balay ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 462a7e14dcfSSatish Balay ierr = VecRestoreArray(Dm,&dm);CHKERRQ(ierr); 463a7e14dcfSSatish Balay } 464a7e14dcfSSatish Balay PetscFunctionReturn(0); 465a7e14dcfSSatish Balay } 466a7e14dcfSSatish Balay 4679fbee547SJacob Faibussowitsch static inline PetscReal ST_InternalPN(PetscScalar in, PetscReal lb, PetscReal ub) 4688370d7cdSHansol Suh { 4698370d7cdSHansol Suh return PetscMax(0,(PetscReal)PetscRealPart(in)-ub) - PetscMax(0,-(PetscReal)PetscRealPart(in)-PetscAbsReal(lb)); 4708370d7cdSHansol Suh } 4718370d7cdSHansol Suh 4729fbee547SJacob Faibussowitsch static inline PetscReal ST_InternalNN(PetscScalar in, PetscReal lb, PetscReal ub) 4738370d7cdSHansol Suh { 4748370d7cdSHansol Suh return PetscMax(0,(PetscReal)PetscRealPart(in) + PetscAbsReal(ub)) - PetscMax(0,-(PetscReal)PetscRealPart(in) - PetscAbsReal(lb)); 4758370d7cdSHansol Suh } 4768370d7cdSHansol Suh 4779fbee547SJacob Faibussowitsch static inline PetscReal ST_InternalPP(PetscScalar in, PetscReal lb, PetscReal ub) 4788370d7cdSHansol Suh { 4798370d7cdSHansol Suh return PetscMax(0, (PetscReal)PetscRealPart(in)-ub) + PetscMin(0, (PetscReal)PetscRealPart(in) - lb); 4808370d7cdSHansol Suh } 4818370d7cdSHansol Suh 4828370d7cdSHansol Suh /*@ 4838370d7cdSHansol Suh TaoSoftThreshold - Calculates soft thresholding routine with input vector 4848370d7cdSHansol Suh and given lower and upper bound and returns it to output vector. 4858370d7cdSHansol Suh 4868370d7cdSHansol Suh Input Parameters: 4878370d7cdSHansol Suh + in - input vector to be thresholded 4888370d7cdSHansol Suh . lb - lower bound 489f0fc11ceSJed Brown - ub - upper bound 4908370d7cdSHansol Suh 49197bb3fdcSJose E. Roman Output Parameter: 4928370d7cdSHansol Suh . out - Soft thresholded output vector 4938370d7cdSHansol Suh 4948370d7cdSHansol Suh Notes: 4958370d7cdSHansol Suh Soft thresholding is defined as 4968370d7cdSHansol Suh \[ S(input,lb,ub) = 4978370d7cdSHansol Suh \begin{cases} 4988370d7cdSHansol Suh input - ub \text{input > ub} \\ 4998370d7cdSHansol Suh 0 \text{lb =< input <= ub} \\ 5008370d7cdSHansol Suh input + lb \text{input < lb} \\ 5018370d7cdSHansol Suh \] 5028370d7cdSHansol Suh 5038370d7cdSHansol Suh Level: developer 5048370d7cdSHansol Suh 5058370d7cdSHansol Suh @*/ 5068370d7cdSHansol Suh PetscErrorCode TaoSoftThreshold(Vec in, PetscReal lb, PetscReal ub, Vec out) 5078370d7cdSHansol Suh { 5088370d7cdSHansol Suh PetscErrorCode ierr; 5098370d7cdSHansol Suh PetscInt i, nlocal, mlocal; 5108370d7cdSHansol Suh PetscScalar *inarray, *outarray; 5118370d7cdSHansol Suh 5128370d7cdSHansol Suh PetscFunctionBegin; 5138370d7cdSHansol Suh ierr = VecGetArrayPair(in, out, &inarray, &outarray);CHKERRQ(ierr); 5148370d7cdSHansol Suh ierr = VecGetLocalSize(in, &nlocal);CHKERRQ(ierr); 5158370d7cdSHansol Suh ierr = VecGetLocalSize(in, &mlocal);CHKERRQ(ierr); 5168370d7cdSHansol Suh 517*3c859ba3SBarry Smith PetscCheck(nlocal == mlocal,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input and output vectors need to be of same size."); 518*3c859ba3SBarry Smith PetscCheck(lb < ub,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Lower bound needs to be lower than upper bound."); 5198370d7cdSHansol Suh 5208370d7cdSHansol Suh if (ub >= 0 && lb < 0) { 5218370d7cdSHansol Suh for (i=0; i<nlocal; i++) outarray[i] = ST_InternalPN(inarray[i], lb, ub); 5228370d7cdSHansol Suh } else if (ub < 0 && lb < 0) { 5238370d7cdSHansol Suh for (i=0; i<nlocal; i++) outarray[i] = ST_InternalNN(inarray[i], lb, ub); 5248370d7cdSHansol Suh } else { 5258370d7cdSHansol Suh for (i=0; i<nlocal; i++) outarray[i] = ST_InternalPP(inarray[i], lb, ub); 5268370d7cdSHansol Suh } 5278370d7cdSHansol Suh 5288370d7cdSHansol Suh ierr = VecRestoreArrayPair(in, out, &inarray, &outarray);CHKERRQ(ierr); 5298370d7cdSHansol Suh PetscFunctionReturn(0); 5308370d7cdSHansol Suh } 531