1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h> 2a7e14dcfSSatish Balay 39371c9d4SSatish Balay PetscErrorCode TaoSetUp_SSFLS(Tao tao) { 4a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay PetscFunctionBegin; 79566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 89566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 99566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->w)); 109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->ff)); 119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->dpsi)); 129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->da)); 139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->db)); 149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->t1)); 159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ssls->t2)); 169566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU)); 17a7e14dcfSSatish Balay PetscFunctionReturn(0); 18a7e14dcfSSatish Balay } 19a7e14dcfSSatish Balay 209371c9d4SSatish Balay static PetscErrorCode TaoSolve_SSFLS(Tao tao) { 21a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 22a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t = 0; 23a7e14dcfSSatish Balay PetscReal delta, rho; 24e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay PetscFunctionBegin; 27a7e14dcfSSatish Balay /* Assume that Setup has been called! 28a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 29a7e14dcfSSatish Balay delta = ssls->delta; 30a7e14dcfSSatish Balay rho = ssls->rho; 31a7e14dcfSSatish Balay 329566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 33a7e14dcfSSatish Balay /* Project solution inside bounds */ 349566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 359566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, Tao_SSLS_FunctionGradient, tao)); 369566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch, Tao_SSLS_Function, tao)); 37a7e14dcfSSatish Balay 38a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 39a7e14dcfSSatish Balay current iterate */ 409566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch, tao->solution, &psi, ssls->dpsi)); 419566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi, NORM_2, &ndpsi)); 42a7e14dcfSSatish Balay 43763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 44d90ca52eSBarry Smith while (PETSC_TRUE) { 4563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "iter: %" PetscInt_FMT ", merit: %g, ndpsi: %g\n", tao->niter, (double)ssls->merit, (double)ndpsi)); 46a7e14dcfSSatish Balay /* Check the termination criteria */ 479566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, ssls->merit, ndpsi, 0.0, tao->ksp_its)); 489566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, ssls->merit, ndpsi, 0.0, t)); 49dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 50763847b4SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) break; 51e1e80dc8SAlp Dener 52e1e80dc8SAlp Dener /* Call general purpose update function */ 53dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 54e6d4cb7fSJason Sarich tao->niter++; 55a7e14dcfSSatish Balay 56a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 57a7e14dcfSSatish Balay rest of the code uses -d.) */ 589566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tao->jacobian, tao->jacobian_pre)); 599566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, ssls->ff, tao->stepdirection)); 609566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &tao->ksp_its)); 61b0026674SJason Sarich tao->ksp_tot_its += tao->ksp_its; 62a7e14dcfSSatish Balay 639566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection, ssls->w)); 649566063dSJacob Faibussowitsch PetscCall(VecScale(ssls->w, -1.0)); 659566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(ssls->w, tao->solution, tao->XL, tao->XU, ssls->w)); 66a7e14dcfSSatish Balay 679566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->w, NORM_2, &normd)); 689566063dSJacob Faibussowitsch PetscCall(VecDot(ssls->w, ssls->dpsi, &innerd)); 69a7e14dcfSSatish Balay 70a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 71d90ca52eSBarry Smith if (innerd >= -delta * PetscPowReal(normd, rho)) { 729566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "newton direction not descent\n")); 739566063dSJacob Faibussowitsch PetscCall(VecCopy(ssls->dpsi, tao->stepdirection)); 749566063dSJacob Faibussowitsch PetscCall(VecDot(ssls->w, ssls->dpsi, &innerd)); 75a7e14dcfSSatish Balay } 76a7e14dcfSSatish Balay 779566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 78a7e14dcfSSatish Balay innerd = -innerd; 79a7e14dcfSSatish Balay 809566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 819566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &psi, ssls->dpsi, tao->stepdirection, &t, &ls_reason)); 829566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi, NORM_2, &ndpsi)); 83a7e14dcfSSatish Balay } 84a7e14dcfSSatish Balay PetscFunctionReturn(0); 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay 879371c9d4SSatish Balay PetscErrorCode TaoDestroy_SSFLS(Tao tao) { 88a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 89a7e14dcfSSatish Balay 90a7e14dcfSSatish Balay PetscFunctionBegin; 919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->ff)); 929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->w)); 939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dpsi)); 949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->da)); 959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->db)); 969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t1)); 979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t2)); 98a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 999566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 100a7e14dcfSSatish Balay PetscFunctionReturn(0); 101a7e14dcfSSatish Balay } 102a7e14dcfSSatish Balay 103a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 1041522df2eSJason Sarich /*MC 1051522df2eSJason Sarich TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving 1061522df2eSJason Sarich complementarity constraints 1071522df2eSJason Sarich 1081522df2eSJason Sarich Options Database Keys: 1091522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 1101522df2eSJason Sarich - -tao_ssls_rho - descent test power 1111522df2eSJason Sarich 1121eb8069cSJason Sarich Level: beginner 1131522df2eSJason Sarich M*/ 1141eb8069cSJason Sarich 1159371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao) { 116a7e14dcfSSatish Balay TAO_SSLS *ssls; 1178caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 118a7e14dcfSSatish Balay 119a7e14dcfSSatish Balay PetscFunctionBegin; 120*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ssls)); 121a7e14dcfSSatish Balay tao->data = (void *)ssls; 122a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_SSFLS; 123a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_SSFLS; 124a7e14dcfSSatish Balay tao->ops->view = TaoView_SSLS; 125a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 126a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSFLS; 127a7e14dcfSSatish Balay 128a7e14dcfSSatish Balay ssls->delta = 1e-10; 129a7e14dcfSSatish Balay ssls->rho = 2.1; 130a7e14dcfSSatish Balay 1319566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 1339566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, armijo_type)); 1349566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix)); 1359566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 13647a47007SBarry Smith /* Linesearch objective and objectivegradient routines are set in solve routine */ 1379566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 1399566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 140a7e14dcfSSatish Balay 1416552cf8aSJason Sarich /* Override default settings (unless already changed) */ 1426552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 1436552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 1446552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 1456552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 1466f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1476552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 1486552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 1496f4723b1SBarry Smith #else 1506552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 1516552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 1526f4723b1SBarry Smith #endif 153a7e14dcfSSatish Balay PetscFunctionReturn(0); 154a7e14dcfSSatish Balay } 155