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; 8*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 9*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 10*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->ff)); 11*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->dpsi)); 12*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->da)); 13*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->db)); 14*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&ssls->t1)); 15*9566063dSJacob 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; 24*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->ff)); 25*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dpsi)); 26*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->da)); 27*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->db)); 28*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t1)); 29*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t2)); 30*9566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 31a7e14dcfSSatish Balay PetscFunctionReturn(0); 32a7e14dcfSSatish Balay } 33a7e14dcfSSatish Balay 34441846f8SBarry Smith static PetscErrorCode TaoSolve_SSILS(Tao tao) 35a7e14dcfSSatish Balay { 36a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 37a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t=0; 38a7e14dcfSSatish Balay PetscReal delta, rho; 39e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay PetscFunctionBegin; 42a7e14dcfSSatish Balay /* Assume that Setup has been called! 43a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 44a7e14dcfSSatish Balay delta = ssls->delta; 45a7e14dcfSSatish Balay rho = ssls->rho; 46a7e14dcfSSatish Balay 47*9566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 48*9566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution)); 49*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao)); 50*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao)); 51a7e14dcfSSatish Balay 52a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 53a7e14dcfSSatish Balay current iterate */ 54*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi)); 55*9566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi)); 56a7e14dcfSSatish Balay 57763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 58d90ca52eSBarry Smith while (PETSC_TRUE) { 59*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi)); 60a7e14dcfSSatish Balay /* Check the termination criteria */ 61*9566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its)); 62*9566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t)); 63*9566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 64763847b4SAlp Dener if (tao->reason!=TAO_CONTINUE_ITERATING) break; 65e1e80dc8SAlp Dener 66e1e80dc8SAlp Dener /* Call general purpose update function */ 67e1e80dc8SAlp Dener if (tao->ops->update) { 68*9566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 69e1e80dc8SAlp Dener } 70e6d4cb7fSJason Sarich tao->niter++; 71a7e14dcfSSatish Balay 72a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 73a7e14dcfSSatish Balay rest of the code uses -d.) */ 74*9566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre)); 75*9566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ssls->ff,tao->stepdirection)); 76*9566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&tao->ksp_its)); 77b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 78*9566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection,NORM_2,&normd)); 79*9566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection,ssls->dpsi,&innerd)); 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 821118d4bcSLisandro Dalcin if (innerd <= delta*PetscPowReal(normd, rho)) { 83*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "newton direction not descent\n")); 84*9566063dSJacob Faibussowitsch PetscCall(VecCopy(ssls->dpsi,tao->stepdirection)); 85*9566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection,ssls->dpsi,&innerd)); 86a7e14dcfSSatish Balay } 87a7e14dcfSSatish Balay 88*9566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 89a7e14dcfSSatish Balay innerd = -innerd; 90a7e14dcfSSatish Balay 91*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 92*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason)); 93*9566063dSJacob Faibussowitsch PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi)); 94a7e14dcfSSatish Balay } 95a7e14dcfSSatish Balay PetscFunctionReturn(0); 96a7e14dcfSSatish Balay } 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 991522df2eSJason Sarich /*MC 1001522df2eSJason Sarich TAOSSILS - semi-smooth infeasible linesearch algorithm for solving 1011522df2eSJason Sarich complementarity constraints 1021522df2eSJason Sarich 1031522df2eSJason Sarich Options Database Keys: 1041522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 1051522df2eSJason Sarich - -tao_ssls_rho - descent test power 1061522df2eSJason Sarich 1071eb8069cSJason Sarich Level: beginner 1081522df2eSJason Sarich M*/ 109728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao) 110a7e14dcfSSatish Balay { 111a7e14dcfSSatish Balay TAO_SSLS *ssls; 1128caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay PetscFunctionBegin; 115*9566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&ssls)); 116a7e14dcfSSatish Balay tao->data = (void*)ssls; 117a7e14dcfSSatish Balay tao->ops->solve=TaoSolve_SSILS; 118a7e14dcfSSatish Balay tao->ops->setup=TaoSetUp_SSILS; 119a7e14dcfSSatish Balay tao->ops->view=TaoView_SSLS; 120a7e14dcfSSatish Balay tao->ops->setfromoptions=TaoSetFromOptions_SSLS; 121a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSILS; 122a7e14dcfSSatish Balay 123a7e14dcfSSatish Balay ssls->delta = 1e-10; 124a7e14dcfSSatish Balay ssls->rho = 2.1; 125a7e14dcfSSatish Balay 126*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch)); 127*9566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 128*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,armijo_type)); 129*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 130*9566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 131a7e14dcfSSatish Balay /* Note: linesearch objective and objectivegradient routines are set in solve routine */ 132*9566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 133*9566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 134*9566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix)); 135a7e14dcfSSatish Balay 1366552cf8aSJason Sarich /* Override default settings (unless already changed) */ 1376552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 1386552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 1396552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 1406552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 1416f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1426552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 1436552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 1446f4723b1SBarry Smith #else 1456552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 1466552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 1476f4723b1SBarry Smith #endif 148a7e14dcfSSatish Balay PetscFunctionReturn(0); 149a7e14dcfSSatish Balay } 150