Actual source code: ssls.c
petsc-3.7.3 2016-08-01
1: #include <../src/tao/complementarity/impls/ssls/ssls.h>
3: /*------------------------------------------------------------*/
6: PetscErrorCode TaoSetFromOptions_SSLS(PetscOptionItems *PetscOptionsObject,Tao tao)
7: {
8: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
12: PetscOptionsHead(PetscOptionsObject,"Semismooth method with a linesearch for complementarity problems");
13: PetscOptionsReal("-ssls_delta", "descent test fraction", "",ssls->delta, &ssls->delta, NULL);
14: PetscOptionsReal("-ssls_rho", "descent test power", "",ssls->rho, &ssls->rho, NULL);
15: TaoLineSearchSetFromOptions(tao->linesearch);
16: KSPSetFromOptions(tao->ksp);
17: PetscOptionsTail();
18: return(0);
19: }
21: /*------------------------------------------------------------*/
24: PetscErrorCode TaoView_SSLS(Tao tao, PetscViewer pv)
25: {
27: return(0);
28: }
30: /*------------------------------------------------------------*/
33: PetscErrorCode Tao_SSLS_Function(TaoLineSearch ls, Vec X, PetscReal *fcn, void *ptr)
34: {
35: Tao tao = (Tao)ptr;
36: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
40: TaoComputeConstraints(tao, X, tao->constraints);
41: VecFischer(X,tao->constraints,tao->XL,tao->XU,ssls->ff);
42: VecNorm(ssls->ff,NORM_2,&ssls->merit);
43: *fcn = 0.5*ssls->merit*ssls->merit;
44: return(0);
45: }
47: /*------------------------------------------------------------*/
50: PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr)
51: {
52: Tao tao = (Tao)ptr;
53: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
57: TaoComputeConstraints(tao, X, tao->constraints);
58: VecFischer(X,tao->constraints,tao->XL,tao->XU,ssls->ff);
59: VecNorm(ssls->ff,NORM_2,&ssls->merit);
60: *fcn = 0.5*ssls->merit*ssls->merit;
62: TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);
64: MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, ssls->t1, ssls->t2,ssls->da, ssls->db);
65: MatDiagonalScale(tao->jacobian,ssls->db,NULL);
66: MatDiagonalSet(tao->jacobian,ssls->da,ADD_VALUES);
67: MatMultTranspose(tao->jacobian,ssls->ff,G);
68: return(0);
69: }