1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h> 2a7e14dcfSSatish Balay 3441846f8SBarry Smith PetscErrorCode TaoSetUp_SSILS(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->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)); 16a7e14dcfSSatish Balay PetscFunctionReturn(0); 17a7e14dcfSSatish Balay } 18a7e14dcfSSatish Balay 19441846f8SBarry Smith PetscErrorCode TaoDestroy_SSILS(Tao tao) 20a7e14dcfSSatish Balay { 21a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->ff)); 259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dpsi)); 269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->da)); 279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->db)); 289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t1)); 299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t2)); 30a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 319566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 32a7e14dcfSSatish Balay PetscFunctionReturn(0); 33a7e14dcfSSatish Balay } 34a7e14dcfSSatish Balay 35441846f8SBarry Smith static PetscErrorCode TaoSolve_SSILS(Tao tao) 36a7e14dcfSSatish Balay { 37a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 38a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t=0; 39a7e14dcfSSatish Balay PetscReal delta, rho; 40e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 41a7e14dcfSSatish Balay 42a7e14dcfSSatish Balay PetscFunctionBegin; 43a7e14dcfSSatish Balay /* Assume that Setup has been called! 44a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 45a7e14dcfSSatish Balay delta = ssls->delta; 46a7e14dcfSSatish Balay rho = ssls->rho; 47a7e14dcfSSatish Balay 489566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 499566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution)); 509566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao)); 519566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao)); 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 54a7e14dcfSSatish Balay current iterate */ 559566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi)); 569566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi)); 57a7e14dcfSSatish Balay 58763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 59d90ca52eSBarry Smith while (PETSC_TRUE) { 6063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "iter: %" PetscInt_FMT ", merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi)); 61a7e14dcfSSatish Balay /* Check the termination criteria */ 629566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its)); 639566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t)); 64*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 65763847b4SAlp Dener if (tao->reason!=TAO_CONTINUE_ITERATING) break; 66e1e80dc8SAlp Dener 67e1e80dc8SAlp Dener /* Call general purpose update function */ 68*dbbe0bcdSBarry Smith PetscTryTypeMethod(tao,update, tao->niter, tao->user_update); 69e6d4cb7fSJason Sarich tao->niter++; 70a7e14dcfSSatish Balay 71a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 72a7e14dcfSSatish Balay rest of the code uses -d.) */ 739566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre)); 749566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ssls->ff,tao->stepdirection)); 759566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&tao->ksp_its)); 76b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 779566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection,NORM_2,&normd)); 789566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection,ssls->dpsi,&innerd)); 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 811118d4bcSLisandro Dalcin if (innerd <= delta*PetscPowReal(normd, rho)) { 829566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "newton direction not descent\n")); 839566063dSJacob Faibussowitsch PetscCall(VecCopy(ssls->dpsi,tao->stepdirection)); 849566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection,ssls->dpsi,&innerd)); 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay 879566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 88a7e14dcfSSatish Balay innerd = -innerd; 89a7e14dcfSSatish Balay 909566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 919566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason)); 929566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi)); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay PetscFunctionReturn(0); 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 981522df2eSJason Sarich /*MC 991522df2eSJason Sarich TAOSSILS - semi-smooth infeasible linesearch algorithm for solving 1001522df2eSJason Sarich complementarity constraints 1011522df2eSJason Sarich 1021522df2eSJason Sarich Options Database Keys: 1031522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 1041522df2eSJason Sarich - -tao_ssls_rho - descent test power 1051522df2eSJason Sarich 1061eb8069cSJason Sarich Level: beginner 1071522df2eSJason Sarich M*/ 108728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao) 109a7e14dcfSSatish Balay { 110a7e14dcfSSatish Balay TAO_SSLS *ssls; 1118caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay PetscFunctionBegin; 1149566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&ssls)); 115a7e14dcfSSatish Balay tao->data = (void*)ssls; 116a7e14dcfSSatish Balay tao->ops->solve=TaoSolve_SSILS; 117a7e14dcfSSatish Balay tao->ops->setup=TaoSetUp_SSILS; 118a7e14dcfSSatish Balay tao->ops->view=TaoView_SSLS; 119a7e14dcfSSatish Balay tao->ops->setfromoptions=TaoSetFromOptions_SSLS; 120a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSILS; 121a7e14dcfSSatish Balay 122a7e14dcfSSatish Balay ssls->delta = 1e-10; 123a7e14dcfSSatish Balay ssls->rho = 2.1; 124a7e14dcfSSatish Balay 1259566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch)); 1269566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 1279566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,armijo_type)); 1289566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 1299566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 130a7e14dcfSSatish Balay /* Note: linesearch objective and objectivegradient routines are set in solve routine */ 1319566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 1339566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix)); 134a7e14dcfSSatish Balay 1356552cf8aSJason Sarich /* Override default settings (unless already changed) */ 1366552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 1376552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 1386552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 1396552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 1406f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1416552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 1426552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 1436f4723b1SBarry Smith #else 1446552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 1456552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 1466f4723b1SBarry Smith #endif 147a7e14dcfSSatish Balay PetscFunctionReturn(0); 148a7e14dcfSSatish Balay } 149