Actual source code: ssls.c

  1: #include <../src/tao/complementarity/impls/ssls/ssls.h>

  3: PetscErrorCode TaoSetFromOptions_SSLS(Tao tao, PetscOptionItems PetscOptionsObject)
  4: {
  5:   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;

  7:   PetscFunctionBegin;
  8:   PetscOptionsHeadBegin(PetscOptionsObject, "Semismooth method with a linesearch for complementarity problems");
  9:   PetscCall(PetscOptionsReal("-ssls_delta", "descent test fraction", "", ssls->delta, &ssls->delta, NULL));
 10:   PetscCall(PetscOptionsReal("-ssls_rho", "descent test power", "", ssls->rho, &ssls->rho, NULL));
 11:   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
 12:   PetscCall(KSPSetFromOptions(tao->ksp));
 13:   PetscOptionsHeadEnd();
 14:   PetscFunctionReturn(PETSC_SUCCESS);
 15: }

 17: PetscErrorCode TaoView_SSLS(Tao tao, PetscViewer pv)
 18: {
 19:   PetscFunctionBegin;
 20:   PetscFunctionReturn(PETSC_SUCCESS);
 21: }

 23: PetscErrorCode Tao_SSLS_Function(TaoLineSearch ls, Vec X, PetscReal *fcn, void *ptr)
 24: {
 25:   Tao       tao  = (Tao)ptr;
 26:   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;

 28:   PetscFunctionBegin;
 29:   PetscCall(TaoComputeConstraints(tao, X, tao->constraints));
 30:   PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, ssls->ff));
 31:   PetscCall(VecNorm(ssls->ff, NORM_2, &ssls->merit));
 32:   *fcn = 0.5 * ssls->merit * ssls->merit;
 33:   PetscFunctionReturn(PETSC_SUCCESS);
 34: }

 36: PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr)
 37: {
 38:   Tao       tao  = (Tao)ptr;
 39:   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;

 41:   PetscFunctionBegin;
 42:   PetscCall(TaoComputeConstraints(tao, X, tao->constraints));
 43:   PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, ssls->ff));
 44:   PetscCall(VecNorm(ssls->ff, NORM_2, &ssls->merit));
 45:   *fcn = 0.5 * ssls->merit * ssls->merit;

 47:   PetscCall(TaoComputeJacobian(tao, tao->solution, tao->jacobian, tao->jacobian_pre));

 49:   PetscCall(MatDFischer(tao->jacobian, tao->solution, tao->constraints, tao->XL, tao->XU, ssls->t1, ssls->t2, ssls->da, ssls->db));
 50:   PetscCall(MatDiagonalScale(tao->jacobian, ssls->db, NULL));
 51:   PetscCall(MatDiagonalSet(tao->jacobian, ssls->da, ADD_VALUES));
 52:   PetscCall(MatMultTranspose(tao->jacobian, ssls->ff, G));
 53:   PetscFunctionReturn(PETSC_SUCCESS);
 54: }