xref: /petsc/src/tao/complementarity/impls/ssls/ssls.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
TaoSetFromOptions_SSLS(Tao tao,PetscOptionItems PetscOptionsObject)3*ce78bad3SBarry Smith PetscErrorCode TaoSetFromOptions_SSLS(Tao tao, PetscOptionItems PetscOptionsObject)
4d71ae5a4SJacob Faibussowitsch {
5a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay   PetscFunctionBegin;
8d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Semismooth method with a linesearch for complementarity problems");
99566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ssls_delta", "descent test fraction", "", ssls->delta, &ssls->delta, NULL));
109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ssls_rho", "descent test power", "", ssls->rho, &ssls->rho, NULL));
119566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
129566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
13d0609cedSBarry Smith   PetscOptionsHeadEnd();
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15a7e14dcfSSatish Balay }
16a7e14dcfSSatish Balay 
TaoView_SSLS(Tao tao,PetscViewer pv)17d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoView_SSLS(Tao tao, PetscViewer pv)
18d71ae5a4SJacob Faibussowitsch {
19a7e14dcfSSatish Balay   PetscFunctionBegin;
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21a7e14dcfSSatish Balay }
22a7e14dcfSSatish Balay 
Tao_SSLS_Function(TaoLineSearch ls,Vec X,PetscReal * fcn,void * ptr)23d71ae5a4SJacob Faibussowitsch PetscErrorCode Tao_SSLS_Function(TaoLineSearch ls, Vec X, PetscReal *fcn, void *ptr)
24d71ae5a4SJacob Faibussowitsch {
25441846f8SBarry Smith   Tao       tao  = (Tao)ptr;
26a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay   PetscFunctionBegin;
299566063dSJacob Faibussowitsch   PetscCall(TaoComputeConstraints(tao, X, tao->constraints));
309566063dSJacob Faibussowitsch   PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, ssls->ff));
319566063dSJacob Faibussowitsch   PetscCall(VecNorm(ssls->ff, NORM_2, &ssls->merit));
32a7e14dcfSSatish Balay   *fcn = 0.5 * ssls->merit * ssls->merit;
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34a7e14dcfSSatish Balay }
35a7e14dcfSSatish Balay 
Tao_SSLS_FunctionGradient(TaoLineSearch ls,Vec X,PetscReal * fcn,Vec G,void * ptr)36d71ae5a4SJacob Faibussowitsch PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr)
37d71ae5a4SJacob Faibussowitsch {
38441846f8SBarry Smith   Tao       tao  = (Tao)ptr;
39a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
40a7e14dcfSSatish Balay 
41a7e14dcfSSatish Balay   PetscFunctionBegin;
429566063dSJacob Faibussowitsch   PetscCall(TaoComputeConstraints(tao, X, tao->constraints));
439566063dSJacob Faibussowitsch   PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, ssls->ff));
449566063dSJacob Faibussowitsch   PetscCall(VecNorm(ssls->ff, NORM_2, &ssls->merit));
45a7e14dcfSSatish Balay   *fcn = 0.5 * ssls->merit * ssls->merit;
46a7e14dcfSSatish Balay 
479566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobian(tao, tao->solution, tao->jacobian, tao->jacobian_pre));
48a7e14dcfSSatish Balay 
499566063dSJacob Faibussowitsch   PetscCall(MatDFischer(tao->jacobian, tao->solution, tao->constraints, tao->XL, tao->XU, ssls->t1, ssls->t2, ssls->da, ssls->db));
509566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(tao->jacobian, ssls->db, NULL));
519566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(tao->jacobian, ssls->da, ADD_VALUES));
529566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian, ssls->ff, G));
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54a7e14dcfSSatish Balay }
55