1 #include <../src/tao/complementarity/impls/ssls/ssls.h> 2 3 /*------------------------------------------------------------*/ 4 PetscErrorCode TaoSetFromOptions_SSLS(Tao tao, PetscOptionItems *PetscOptionsObject) 5 { 6 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 7 8 PetscFunctionBegin; 9 PetscOptionsHeadBegin(PetscOptionsObject, "Semismooth method with a linesearch for complementarity problems"); 10 PetscCall(PetscOptionsReal("-ssls_delta", "descent test fraction", "", ssls->delta, &ssls->delta, NULL)); 11 PetscCall(PetscOptionsReal("-ssls_rho", "descent test power", "", ssls->rho, &ssls->rho, NULL)); 12 PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 13 PetscCall(KSPSetFromOptions(tao->ksp)); 14 PetscOptionsHeadEnd(); 15 PetscFunctionReturn(0); 16 } 17 18 /*------------------------------------------------------------*/ 19 PetscErrorCode TaoView_SSLS(Tao tao, PetscViewer pv) 20 { 21 PetscFunctionBegin; 22 PetscFunctionReturn(0); 23 } 24 25 /*------------------------------------------------------------*/ 26 PetscErrorCode Tao_SSLS_Function(TaoLineSearch ls, Vec X, PetscReal *fcn, void *ptr) 27 { 28 Tao tao = (Tao)ptr; 29 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 30 31 PetscFunctionBegin; 32 PetscCall(TaoComputeConstraints(tao, X, tao->constraints)); 33 PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, ssls->ff)); 34 PetscCall(VecNorm(ssls->ff, NORM_2, &ssls->merit)); 35 *fcn = 0.5 * ssls->merit * ssls->merit; 36 PetscFunctionReturn(0); 37 } 38 39 /*------------------------------------------------------------*/ 40 PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr) 41 { 42 Tao tao = (Tao)ptr; 43 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 44 45 PetscFunctionBegin; 46 PetscCall(TaoComputeConstraints(tao, X, tao->constraints)); 47 PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, ssls->ff)); 48 PetscCall(VecNorm(ssls->ff, NORM_2, &ssls->merit)); 49 *fcn = 0.5 * ssls->merit * ssls->merit; 50 51 PetscCall(TaoComputeJacobian(tao, tao->solution, tao->jacobian, tao->jacobian_pre)); 52 53 PetscCall(MatDFischer(tao->jacobian, tao->solution, tao->constraints, tao->XL, tao->XU, ssls->t1, ssls->t2, ssls->da, ssls->db)); 54 PetscCall(MatDiagonalScale(tao->jacobian, ssls->db, NULL)); 55 PetscCall(MatDiagonalSet(tao->jacobian, ssls->da, ADD_VALUES)); 56 PetscCall(MatMultTranspose(tao->jacobian, ssls->ff, G)); 57 PetscFunctionReturn(0); 58 } 59