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