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