1a7e14dcfSSatish Balay /* Context for SSXLS 2a7e14dcfSSatish Balay -- semismooth (SS) - function not differentiable 3a7e14dcfSSatish Balay - merit function continuously differentiable 4a7e14dcfSSatish Balay - Fischer-Burmeister reformulation of complementarity 5a7e14dcfSSatish Balay - Billups composition for two finite bounds 6a7e14dcfSSatish Balay -- infeasible (I) - iterates not guaranteed to remain within bounds 7a7e14dcfSSatish Balay -- feasible (F) - iterates guaranteed to remain within bounds 8a7e14dcfSSatish Balay -- linesearch (LS) - Armijo rule on direction 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay Many other reformulations are possible and combinations of 11a7e14dcfSSatish Balay feasible/infeasible and linesearch/trust region are possible. 12a7e14dcfSSatish Balay 13a7e14dcfSSatish Balay Basic theory 14a7e14dcfSSatish Balay Fischer-Burmeister reformulation is semismooth with a continuously 15a7e14dcfSSatish Balay differentiable merit function and strongly semismooth if the F has 16a7e14dcfSSatish Balay lipschitz continuous derivatives. 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay Every accumulation point generated by the algorithm is a stationary 19a7e14dcfSSatish Balay point for the merit function. Stationary points of the merit function 20a7e14dcfSSatish Balay are solutions of the complementarity problem if 21a7e14dcfSSatish Balay a. the stationary point has a BD-regular subdifferential, or 22a7e14dcfSSatish Balay b. the Schur complement F'/F'_ff is a P_0-matrix where ff is the 23a7e14dcfSSatish Balay index set corresponding to the free variables. 24a7e14dcfSSatish Balay 25a7e14dcfSSatish Balay If one of the accumulation points has a BD-regular subdifferential then 26a7e14dcfSSatish Balay a. the entire sequence converges to this accumulation point at 27a7e14dcfSSatish Balay a local q-superlinear rate 28a7e14dcfSSatish Balay b. if in addition the reformulation is strongly semismooth near 29a7e14dcfSSatish Balay this accumulation point, then the algorithm converges at a 30a7e14dcfSSatish Balay local q-quadratic rate. 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay The theory for the feasible version follows from the feasible descent 33a7e14dcfSSatish Balay algorithm framework. 34a7e14dcfSSatish Balay 35a7e14dcfSSatish Balay References: 36606c0280SSatish Balay + * - Billups, "Algorithms for Complementarity Problems and Generalized 37a7e14dcfSSatish Balay Equations," Ph.D thesis, University of Wisconsin - Madison, 1995. 38606c0280SSatish Balay . * - De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the 39a7e14dcfSSatish Balay Solution of Nonlinear Complementarity Problems," Mathematical 40a7e14dcfSSatish Balay Programming, 75, pages 407-439, 1996. 41606c0280SSatish Balay . * - Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed 42a7e14dcfSSatish Balay Complementarity Problems," Mathematical Programming, 86, 43a7e14dcfSSatish Balay pages 475-497, 1999. 44606c0280SSatish Balay . * - Fischer, "A Special Newton-type Optimization Method," Optimization, 45a7e14dcfSSatish Balay 24, pages 269-284, 1992 46606c0280SSatish Balay - * - Munson, Facchinei, Ferris, Fischer, Kanzow, "The Semismooth Algorithm 47a7e14dcfSSatish Balay for Large Scale Complementarity Problems," Technical Report 99-06, 48a7e14dcfSSatish Balay University of Wisconsin - Madison, 1999. 49a7e14dcfSSatish Balay */ 50a7e14dcfSSatish Balay 51a7e14dcfSSatish Balay #ifndef __TAO_SSLS_H 52a7e14dcfSSatish Balay #define __TAO_SSLS_H 53af0996ceSBarry Smith #include <petsc/private/taoimpl.h> 54a7e14dcfSSatish Balay 55a7e14dcfSSatish Balay typedef struct { 56a7e14dcfSSatish Balay Vec ff; /* fischer function */ 57a7e14dcfSSatish Balay Vec dpsi; /* gradient of psi */ 58a7e14dcfSSatish Balay 59a7e14dcfSSatish Balay Vec da; /* work vector for subdifferential calculation (diag pert) */ 60a7e14dcfSSatish Balay Vec db; /* work vector for subdifferential calculation (row scale) */ 61a7e14dcfSSatish Balay Vec dm; /* work vector for subdifferential calculation (mu vector) */ 62a7e14dcfSSatish Balay Vec dxfree; 63a7e14dcfSSatish Balay 64a7e14dcfSSatish Balay Vec t1; /* work vector */ 65a7e14dcfSSatish Balay Vec t2; /* work vector */ 66a7e14dcfSSatish Balay 67a7e14dcfSSatish Balay Vec r1,r2,r3,w; /* work vectors */ 68a7e14dcfSSatish Balay 69a7e14dcfSSatish Balay PetscReal merit; /* merit function value (norm(fischer)) */ 70a7e14dcfSSatish Balay PetscReal merit_eqn; 71a7e14dcfSSatish Balay PetscReal merit_mu; 72a7e14dcfSSatish Balay 73a7e14dcfSSatish Balay PetscReal delta; 74a7e14dcfSSatish Balay PetscReal rho; 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay PetscReal rtol; /* Solution tolerances */ 77a7e14dcfSSatish Balay PetscReal atol; 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay PetscReal identifier; /* Active-set identification */ 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay /* Interior-point method data */ 82a7e14dcfSSatish Balay PetscReal mu_init; /* initial smoothing parameter value */ 83a7e14dcfSSatish Balay PetscReal mu; /* smoothing parameter */ 84a7e14dcfSSatish Balay PetscReal dmu; /* direction in smoothing parameter */ 85a7e14dcfSSatish Balay PetscReal mucon; /* smoothing parameter constraint */ 86a7e14dcfSSatish Balay PetscReal d_mucon; /* derivative of smoothing constraint with respect to mu */ 87a7e14dcfSSatish Balay PetscReal g_mucon; /* gradient of merit function with respect to mu */ 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay Mat J_sub, Jpre_sub; /* subset of jacobian */ 90a7e14dcfSSatish Balay Vec f; /* constraint function */ 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay IS fixed; 93a7e14dcfSSatish Balay IS free; 94a7e14dcfSSatish Balay } TAO_SSLS; 95a7e14dcfSSatish Balay 96*dbbe0bcdSBarry Smith PetscErrorCode TaoSetFromOptions_SSLS(Tao,PetscOptionItems *); 97441846f8SBarry Smith PetscErrorCode TaoView_SSLS(Tao,PetscViewer); 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay PetscErrorCode Tao_SSLS_Function(TaoLineSearch, Vec, PetscReal *, void *); 100a7e14dcfSSatish Balay PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *); 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay #endif 103a7e14dcfSSatish Balay 104