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