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 331d27aa22SBarry Smith algorithm framework. See {cite}`billups:algorithms`, {cite}`deluca.facchinei.ea:semismooth`, 341d27aa22SBarry Smith {cite}`ferris.kanzow.ea:feasible`, {cite}`fischer:special`, and {cite}`munson.facchinei.ea:semismooth`. 35a7e14dcfSSatish Balay */ 36a7e14dcfSSatish Balay 37a4963045SJacob Faibussowitsch #pragma once 38af0996ceSBarry Smith #include <petsc/private/taoimpl.h> 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay typedef struct { 41a7e14dcfSSatish Balay Vec ff; /* fischer function */ 42a7e14dcfSSatish Balay Vec dpsi; /* gradient of psi */ 43a7e14dcfSSatish Balay 44a7e14dcfSSatish Balay Vec da; /* work vector for subdifferential calculation (diag pert) */ 45a7e14dcfSSatish Balay Vec db; /* work vector for subdifferential calculation (row scale) */ 46a7e14dcfSSatish Balay Vec dm; /* work vector for subdifferential calculation (mu vector) */ 47a7e14dcfSSatish Balay Vec dxfree; 48a7e14dcfSSatish Balay 49a7e14dcfSSatish Balay Vec t1; /* work vector */ 50a7e14dcfSSatish Balay Vec t2; /* work vector */ 51a7e14dcfSSatish Balay 52a7e14dcfSSatish Balay Vec r1, r2, r3, w; /* work vectors */ 53a7e14dcfSSatish Balay 54a7e14dcfSSatish Balay PetscReal merit; /* merit function value (norm(fischer)) */ 55a7e14dcfSSatish Balay PetscReal merit_eqn; 56a7e14dcfSSatish Balay PetscReal merit_mu; 57a7e14dcfSSatish Balay 58a7e14dcfSSatish Balay PetscReal delta; 59a7e14dcfSSatish Balay PetscReal rho; 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay PetscReal rtol; /* Solution tolerances */ 62a7e14dcfSSatish Balay PetscReal atol; 63a7e14dcfSSatish Balay 64a7e14dcfSSatish Balay PetscReal identifier; /* Active-set identification */ 65a7e14dcfSSatish Balay 66a7e14dcfSSatish Balay /* Interior-point method data */ 67a7e14dcfSSatish Balay PetscReal mu_init; /* initial smoothing parameter value */ 68a7e14dcfSSatish Balay PetscReal mu; /* smoothing parameter */ 69a7e14dcfSSatish Balay PetscReal dmu; /* direction in smoothing parameter */ 70a7e14dcfSSatish Balay PetscReal mucon; /* smoothing parameter constraint */ 71a7e14dcfSSatish Balay PetscReal d_mucon; /* derivative of smoothing constraint with respect to mu */ 72a7e14dcfSSatish Balay PetscReal g_mucon; /* gradient of merit function with respect to mu */ 73a7e14dcfSSatish Balay 74a7e14dcfSSatish Balay Mat J_sub, Jpre_sub; /* subset of jacobian */ 75a7e14dcfSSatish Balay Vec f; /* constraint function */ 76a7e14dcfSSatish Balay 77a7e14dcfSSatish Balay IS fixed; 78a7e14dcfSSatish Balay IS free; 79a7e14dcfSSatish Balay } TAO_SSLS; 80a7e14dcfSSatish Balay 81*ce78bad3SBarry Smith PetscErrorCode TaoSetFromOptions_SSLS(Tao, PetscOptionItems); 82441846f8SBarry Smith PetscErrorCode TaoView_SSLS(Tao, PetscViewer); 83a7e14dcfSSatish Balay 84a7e14dcfSSatish Balay PetscErrorCode Tao_SSLS_Function(TaoLineSearch, Vec, PetscReal *, void *); 85a7e14dcfSSatish Balay PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *); 86