1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h> 2a7e14dcfSSatish Balay 3a7e14dcfSSatish Balay #undef __FUNCT__ 4a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_SSFLS" 5441846f8SBarry Smith PetscErrorCode TaoSetUp_SSFLS(Tao tao) 6a7e14dcfSSatish Balay { 7a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 8a7e14dcfSSatish Balay PetscErrorCode ierr; 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay PetscFunctionBegin; 11a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 12a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 13a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->w);CHKERRQ(ierr); 14a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->ff);CHKERRQ(ierr); 15a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr); 16a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr); 17a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr); 18a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr); 19a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr); 20a7e14dcfSSatish Balay if (!tao->XL) { 21a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr); 22a7e14dcfSSatish Balay } 23a7e14dcfSSatish Balay if (!tao->XU) { 24a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr); 25a7e14dcfSSatish Balay } 26a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 27a7e14dcfSSatish Balay PetscFunctionReturn(0); 28a7e14dcfSSatish Balay } 29a7e14dcfSSatish Balay 30a7e14dcfSSatish Balay #undef __FUNCT__ 31a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_SSFLS" 32441846f8SBarry Smith static PetscErrorCode TaoSolve_SSFLS(Tao tao) 33a7e14dcfSSatish Balay { 34a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 35a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t=0; 36a7e14dcfSSatish Balay PetscReal delta, rho; 37e4cb33bbSBarry Smith TaoConvergedReason reason; 38e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 39a7e14dcfSSatish Balay PetscErrorCode ierr; 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 47a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 48a7e14dcfSSatish Balay /* Project solution inside bounds */ 49a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 5047a47007SBarry Smith ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr); 51a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 54a7e14dcfSSatish Balay current iterate */ 55a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr); 56a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 57a7e14dcfSSatish Balay 58d90ca52eSBarry Smith while (PETSC_TRUE) { 59*8931d482SJason Sarich ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr); 60a7e14dcfSSatish Balay /* Check the termination criteria */ 61*8931d482SJason Sarich ierr = TaoMonitor(tao,tao->niter++,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr); 62a7e14dcfSSatish Balay if (reason!=TAO_CONTINUE_ITERATING) break; 63a7e14dcfSSatish Balay 64a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 65a7e14dcfSSatish Balay rest of the code uses -d.) */ 6623ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr); 67a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr); 68b0026674SJason Sarich ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr); 69b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 70a7e14dcfSSatish Balay 71a7e14dcfSSatish Balay ierr = VecCopy(tao->stepdirection,ssls->w);CHKERRQ(ierr); 72a7e14dcfSSatish Balay ierr = VecScale(ssls->w,-1.0);CHKERRQ(ierr); 73a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);CHKERRQ(ierr); 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay ierr = VecNorm(ssls->w,NORM_2,&normd);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 77a7e14dcfSSatish Balay 78a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 79d90ca52eSBarry Smith if (innerd >= -delta*PetscPowReal(normd, rho)) { 80a7e14dcfSSatish Balay ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr); 81a7e14dcfSSatish Balay ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr); 82a7e14dcfSSatish Balay ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 83a7e14dcfSSatish Balay } 84a7e14dcfSSatish Balay 85a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 86a7e14dcfSSatish Balay innerd = -innerd; 87a7e14dcfSSatish Balay 88302440fdSBarry Smith ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 89a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr); 90a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 91a7e14dcfSSatish Balay } 92a7e14dcfSSatish Balay PetscFunctionReturn(0); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay 95a7e14dcfSSatish Balay #undef __FUNCT__ 96a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_SSFLS" 97441846f8SBarry Smith PetscErrorCode TaoDestroy_SSFLS(Tao tao) 98a7e14dcfSSatish Balay { 99a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 100a7e14dcfSSatish Balay PetscErrorCode ierr; 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay PetscFunctionBegin; 103a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr); 104a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->w);CHKERRQ(ierr); 105a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr); 106a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->da);CHKERRQ(ierr); 107a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->db);CHKERRQ(ierr); 108a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr); 109a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr); 110a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 111a7e14dcfSSatish Balay PetscFunctionReturn(0); 112a7e14dcfSSatish Balay } 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 1151522df2eSJason Sarich /*MC 1161522df2eSJason Sarich TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving 1171522df2eSJason Sarich complementarity constraints 1181522df2eSJason Sarich 1191522df2eSJason Sarich Options Database Keys: 1201522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 1211522df2eSJason Sarich - -tao_ssls_rho - descent test power 1221522df2eSJason Sarich 1231eb8069cSJason Sarich Level: beginner 1241522df2eSJason Sarich M*/ 1251eb8069cSJason Sarich 126a7e14dcfSSatish Balay #undef __FUNCT__ 127a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_SSFLS" 128728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao) 129a7e14dcfSSatish Balay { 130a7e14dcfSSatish Balay TAO_SSLS *ssls; 131a7e14dcfSSatish Balay PetscErrorCode ierr; 1328caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay PetscFunctionBegin; 1353c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr); 136a7e14dcfSSatish Balay tao->data = (void*)ssls; 137a7e14dcfSSatish Balay tao->ops->solve=TaoSolve_SSFLS; 138a7e14dcfSSatish Balay tao->ops->setup=TaoSetUp_SSFLS; 139a7e14dcfSSatish Balay tao->ops->view=TaoView_SSLS; 140a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 141a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSFLS; 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay ssls->delta = 1e-10; 144a7e14dcfSSatish Balay ssls->rho = 2.1; 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 147a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr); 148302440fdSBarry Smith ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 14947a47007SBarry Smith /* Linesearch objective and objectivegradient routines are set in solve routine */ 150a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 151a7e14dcfSSatish Balay 152a7e14dcfSSatish Balay tao->max_it = 2000; 153a7e14dcfSSatish Balay tao->max_funcs = 4000; 1546f4723b1SBarry Smith tao->fatol = 0; 1556f4723b1SBarry Smith tao->frtol = 0; 1566f4723b1SBarry Smith tao->grtol=0; 1576f4723b1SBarry Smith tao->grtol=0; 1586f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1596f4723b1SBarry Smith tao->gatol = 1.0e-6; 1606f4723b1SBarry Smith tao->fmin = 1.0e-4; 1616f4723b1SBarry Smith #else 162a7e14dcfSSatish Balay tao->gatol = 1.0e-16; 163a7e14dcfSSatish Balay tao->fmin = 1.0e-8; 1646f4723b1SBarry Smith #endif 165a7e14dcfSSatish Balay PetscFunctionReturn(0); 166a7e14dcfSSatish Balay } 167728e0ed0SBarry Smith 168a7e14dcfSSatish Balay 169