1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h> 2a7e14dcfSSatish Balay 3441846f8SBarry Smith PetscErrorCode TaoSetUp_SSFLS(Tao tao) 4a7e14dcfSSatish Balay { 5a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 6a7e14dcfSSatish Balay 7a7e14dcfSSatish Balay PetscFunctionBegin; 89566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 99566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->w)); 119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->ff)); 129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->dpsi)); 139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->da)); 149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->db)); 159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->t1)); 169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->t2)); 179566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU)); 18a7e14dcfSSatish Balay PetscFunctionReturn(0); 19a7e14dcfSSatish Balay } 20a7e14dcfSSatish Balay 21441846f8SBarry Smith static PetscErrorCode TaoSolve_SSFLS(Tao tao) 22a7e14dcfSSatish Balay { 23a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 24a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t=0; 25a7e14dcfSSatish Balay PetscReal delta, rho; 26e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 27a7e14dcfSSatish Balay 28a7e14dcfSSatish Balay PetscFunctionBegin; 29a7e14dcfSSatish Balay /* Assume that Setup has been called! 30a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 31a7e14dcfSSatish Balay delta = ssls->delta; 32a7e14dcfSSatish Balay rho = ssls->rho; 33a7e14dcfSSatish Balay 349566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 35a7e14dcfSSatish Balay /* Project solution inside bounds */ 369566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution)); 379566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao)); 389566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao)); 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 41a7e14dcfSSatish Balay current iterate */ 429566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi)); 439566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi)); 44a7e14dcfSSatish Balay 45763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 46d90ca52eSBarry Smith while (PETSC_TRUE) { 4763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "iter: %" PetscInt_FMT ", merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi)); 48a7e14dcfSSatish Balay /* Check the termination criteria */ 499566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its)); 509566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t)); 519566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 52763847b4SAlp Dener if (tao->reason!=TAO_CONTINUE_ITERATING) break; 53e1e80dc8SAlp Dener 54e1e80dc8SAlp Dener /* Call general purpose update function */ 55e1e80dc8SAlp Dener if (tao->ops->update) { 569566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 57e1e80dc8SAlp Dener } 58e6d4cb7fSJason Sarich tao->niter++; 59a7e14dcfSSatish Balay 60a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 61a7e14dcfSSatish Balay rest of the code uses -d.) */ 629566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre)); 639566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ssls->ff,tao->stepdirection)); 649566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&tao->ksp_its)); 65b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 66a7e14dcfSSatish Balay 679566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection,ssls->w)); 689566063dSJacob Faibussowitsch PetscCall(VecScale(ssls->w,-1.0)); 699566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w)); 70a7e14dcfSSatish Balay 719566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->w,NORM_2,&normd)); 729566063dSJacob Faibussowitsch PetscCall(VecDot(ssls->w,ssls->dpsi,&innerd)); 73a7e14dcfSSatish Balay 74a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 75d90ca52eSBarry Smith if (innerd >= -delta*PetscPowReal(normd, rho)) { 769566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "newton direction not descent\n")); 779566063dSJacob Faibussowitsch PetscCall(VecCopy(ssls->dpsi,tao->stepdirection)); 789566063dSJacob Faibussowitsch PetscCall(VecDot(ssls->w,ssls->dpsi,&innerd)); 79a7e14dcfSSatish Balay } 80a7e14dcfSSatish Balay 819566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 82a7e14dcfSSatish Balay innerd = -innerd; 83a7e14dcfSSatish Balay 849566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 859566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason)); 869566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi)); 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay PetscFunctionReturn(0); 89a7e14dcfSSatish Balay } 90a7e14dcfSSatish Balay 91441846f8SBarry Smith PetscErrorCode TaoDestroy_SSFLS(Tao tao) 92a7e14dcfSSatish Balay { 93a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 94a7e14dcfSSatish Balay 95a7e14dcfSSatish Balay PetscFunctionBegin; 969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->ff)); 979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->w)); 989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dpsi)); 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->da)); 1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->db)); 1019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t1)); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t2)); 103*a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 105a7e14dcfSSatish Balay PetscFunctionReturn(0); 106a7e14dcfSSatish Balay } 107a7e14dcfSSatish Balay 108a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 1091522df2eSJason Sarich /*MC 1101522df2eSJason Sarich TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving 1111522df2eSJason Sarich complementarity constraints 1121522df2eSJason Sarich 1131522df2eSJason Sarich Options Database Keys: 1141522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 1151522df2eSJason Sarich - -tao_ssls_rho - descent test power 1161522df2eSJason Sarich 1171eb8069cSJason Sarich Level: beginner 1181522df2eSJason Sarich M*/ 1191eb8069cSJason Sarich 120728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao) 121a7e14dcfSSatish Balay { 122a7e14dcfSSatish Balay TAO_SSLS *ssls; 1238caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay PetscFunctionBegin; 1269566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&ssls)); 127a7e14dcfSSatish Balay tao->data = (void*)ssls; 128a7e14dcfSSatish Balay tao->ops->solve=TaoSolve_SSFLS; 129a7e14dcfSSatish Balay tao->ops->setup=TaoSetUp_SSFLS; 130a7e14dcfSSatish Balay tao->ops->view=TaoView_SSLS; 131a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 132a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSFLS; 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay ssls->delta = 1e-10; 135a7e14dcfSSatish Balay ssls->rho = 2.1; 136a7e14dcfSSatish Balay 1379566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch)); 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 1399566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,armijo_type)); 1409566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 1419566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 14247a47007SBarry Smith /* Linesearch objective and objectivegradient routines are set in solve routine */ 1439566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 1449566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 1459566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix)); 146a7e14dcfSSatish Balay 1476552cf8aSJason Sarich /* Override default settings (unless already changed) */ 1486552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 1496552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 1506552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 1516552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 1526f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1536552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 1546552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 1556f4723b1SBarry Smith #else 1566552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 1576552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 1586f4723b1SBarry Smith #endif 159a7e14dcfSSatish Balay PetscFunctionReturn(0); 160a7e14dcfSSatish Balay } 161