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