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 PetscErrorCode ierr; 7a7e14dcfSSatish Balay 8a7e14dcfSSatish Balay PetscFunctionBegin; 9a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 10a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 11a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->w);CHKERRQ(ierr); 12a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->ff);CHKERRQ(ierr); 13a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr); 14a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr); 15a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr); 16a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr); 17a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr); 18a7e14dcfSSatish Balay if (!tao->XL) { 19a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr); 20a7e14dcfSSatish Balay } 21a7e14dcfSSatish Balay if (!tao->XU) { 22a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr); 23a7e14dcfSSatish Balay } 24a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 25a7e14dcfSSatish Balay PetscFunctionReturn(0); 26a7e14dcfSSatish Balay } 27a7e14dcfSSatish Balay 28441846f8SBarry Smith static PetscErrorCode TaoSolve_SSFLS(Tao tao) 29a7e14dcfSSatish Balay { 30a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 31a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t=0; 32a7e14dcfSSatish Balay PetscReal delta, rho; 33e4cb33bbSBarry Smith TaoConvergedReason reason; 34e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 35a7e14dcfSSatish Balay PetscErrorCode ierr; 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay PetscFunctionBegin; 38a7e14dcfSSatish Balay /* Assume that Setup has been called! 39a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 40a7e14dcfSSatish Balay delta = ssls->delta; 41a7e14dcfSSatish Balay rho = ssls->rho; 42a7e14dcfSSatish Balay 43a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 44a7e14dcfSSatish Balay /* Project solution inside bounds */ 45a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 4647a47007SBarry Smith ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr); 47a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 48a7e14dcfSSatish Balay 49a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 50a7e14dcfSSatish Balay current iterate */ 51a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr); 52a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 53a7e14dcfSSatish Balay 54d90ca52eSBarry Smith while (PETSC_TRUE) { 558931d482SJason Sarich ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr); 56a7e14dcfSSatish Balay /* Check the termination criteria */ 57e6d4cb7fSJason Sarich ierr = TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr); 58a7e14dcfSSatish Balay if (reason!=TAO_CONTINUE_ITERATING) break; 59e6d4cb7fSJason Sarich tao->niter++; 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 62a7e14dcfSSatish Balay rest of the code uses -d.) */ 6323ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr); 65b0026674SJason Sarich ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr); 66b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay ierr = VecCopy(tao->stepdirection,ssls->w);CHKERRQ(ierr); 69a7e14dcfSSatish Balay ierr = VecScale(ssls->w,-1.0);CHKERRQ(ierr); 70a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);CHKERRQ(ierr); 71a7e14dcfSSatish Balay 72a7e14dcfSSatish Balay ierr = VecNorm(ssls->w,NORM_2,&normd);CHKERRQ(ierr); 73a7e14dcfSSatish Balay ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 76d90ca52eSBarry Smith if (innerd >= -delta*PetscPowReal(normd, rho)) { 77a7e14dcfSSatish Balay ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr); 79a7e14dcfSSatish Balay ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 80a7e14dcfSSatish Balay } 81a7e14dcfSSatish Balay 82a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 83a7e14dcfSSatish Balay innerd = -innerd; 84a7e14dcfSSatish Balay 85302440fdSBarry Smith ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 86a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 88a7e14dcfSSatish Balay } 89a7e14dcfSSatish Balay PetscFunctionReturn(0); 90a7e14dcfSSatish Balay } 91a7e14dcfSSatish Balay 92441846f8SBarry Smith PetscErrorCode TaoDestroy_SSFLS(Tao tao) 93a7e14dcfSSatish Balay { 94a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 95a7e14dcfSSatish Balay PetscErrorCode ierr; 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay PetscFunctionBegin; 98a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr); 99a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->w);CHKERRQ(ierr); 100a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr); 101a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->da);CHKERRQ(ierr); 102a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->db);CHKERRQ(ierr); 103a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr); 104a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr); 105a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 106a7e14dcfSSatish Balay PetscFunctionReturn(0); 107a7e14dcfSSatish Balay } 108a7e14dcfSSatish Balay 109a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 1101522df2eSJason Sarich /*MC 1111522df2eSJason Sarich TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving 1121522df2eSJason Sarich complementarity constraints 1131522df2eSJason Sarich 1141522df2eSJason Sarich Options Database Keys: 1151522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 1161522df2eSJason Sarich - -tao_ssls_rho - descent test power 1171522df2eSJason Sarich 1181eb8069cSJason Sarich Level: beginner 1191522df2eSJason Sarich M*/ 1201eb8069cSJason Sarich 121728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao) 122a7e14dcfSSatish Balay { 123a7e14dcfSSatish Balay TAO_SSLS *ssls; 124a7e14dcfSSatish Balay PetscErrorCode ierr; 1258caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 126a7e14dcfSSatish Balay 127a7e14dcfSSatish Balay PetscFunctionBegin; 1283c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr); 129a7e14dcfSSatish Balay tao->data = (void*)ssls; 130a7e14dcfSSatish Balay tao->ops->solve=TaoSolve_SSFLS; 131a7e14dcfSSatish Balay tao->ops->setup=TaoSetUp_SSFLS; 132a7e14dcfSSatish Balay tao->ops->view=TaoView_SSLS; 133a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 134a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSFLS; 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay ssls->delta = 1e-10; 137a7e14dcfSSatish Balay ssls->rho = 2.1; 138a7e14dcfSSatish Balay 139a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 140*63b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 141a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr); 1425d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 143302440fdSBarry Smith ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 14447a47007SBarry Smith /* Linesearch objective and objectivegradient routines are set in solve routine */ 145a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 146*63b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1475d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 148a7e14dcfSSatish Balay 1496552cf8aSJason Sarich /* Override default settings (unless already changed) */ 1506552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 1516552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 1526552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 1536552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 1546f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1556552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 1566552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 1576f4723b1SBarry Smith #else 1586552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 1596552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 1606f4723b1SBarry Smith #endif 161a7e14dcfSSatish Balay PetscFunctionReturn(0); 162a7e14dcfSSatish Balay } 163