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