Actual source code: ssls.c
petsc-3.9.4 2018-09-11
1: #include <../src/tao/complementarity/impls/ssls/ssls.h>
3: /*------------------------------------------------------------*/
4: PetscErrorCode TaoSetFromOptions_SSLS(PetscOptionItems *PetscOptionsObject,Tao tao)
5: {
6: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
10: PetscOptionsHead(PetscOptionsObject,"Semismooth method with a linesearch for complementarity problems");
11: PetscOptionsReal("-ssls_delta", "descent test fraction", "",ssls->delta, &ssls->delta, NULL);
12: PetscOptionsReal("-ssls_rho", "descent test power", "",ssls->rho, &ssls->rho, NULL);
13: TaoLineSearchSetFromOptions(tao->linesearch);
14: KSPSetFromOptions(tao->ksp);
15: PetscOptionsTail();
16: return(0);
17: }
19: /*------------------------------------------------------------*/
20: PetscErrorCode TaoView_SSLS(Tao tao, PetscViewer pv)
21: {
23: return(0);
24: }
26: /*------------------------------------------------------------*/
27: PetscErrorCode Tao_SSLS_Function(TaoLineSearch ls, Vec X, PetscReal *fcn, void *ptr)
28: {
29: Tao tao = (Tao)ptr;
30: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
34: TaoComputeConstraints(tao, X, tao->constraints);
35: VecFischer(X,tao->constraints,tao->XL,tao->XU,ssls->ff);
36: VecNorm(ssls->ff,NORM_2,&ssls->merit);
37: *fcn = 0.5*ssls->merit*ssls->merit;
38: return(0);
39: }
41: /*------------------------------------------------------------*/
42: PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr)
43: {
44: Tao tao = (Tao)ptr;
45: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
49: TaoComputeConstraints(tao, X, tao->constraints);
50: VecFischer(X,tao->constraints,tao->XL,tao->XU,ssls->ff);
51: VecNorm(ssls->ff,NORM_2,&ssls->merit);
52: *fcn = 0.5*ssls->merit*ssls->merit;
54: TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);
56: MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, ssls->t1, ssls->t2,ssls->da, ssls->db);
57: MatDiagonalScale(tao->jacobian,ssls->db,NULL);
58: MatDiagonalSet(tao->jacobian,ssls->da,ADD_VALUES);
59: MatMultTranspose(tao->jacobian,ssls->ff,G);
60: return(0);
61: }