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; 37a7e14dcfSSatish Balay PetscInt iter=0,kspits; 38*e4cb33bbSBarry Smith TaoConvergedReason reason; 39*e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 40a7e14dcfSSatish Balay PetscErrorCode ierr; 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 48a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 49a7e14dcfSSatish Balay /* Project solution inside bounds */ 50a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 5147a47007SBarry Smith ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr); 52a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 53a7e14dcfSSatish Balay 54a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 55a7e14dcfSSatish Balay current iterate */ 56a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr); 57a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 58a7e14dcfSSatish Balay 59a7e14dcfSSatish Balay while (1) { 6047a47007SBarry Smith ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",iter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr); 61a7e14dcfSSatish Balay /* Check the termination criteria */ 6247a47007SBarry Smith ierr = TaoMonitor(tao,iter++,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr); 63a7e14dcfSSatish Balay if (reason!=TAO_CONTINUE_ITERATING) break; 64a7e14dcfSSatish Balay 65a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 66a7e14dcfSSatish Balay rest of the code uses -d.) */ 67a7e14dcfSSatish Balay ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre,ssls->matflag);CHKERRQ(ierr); 68a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr); 69a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 70a7e14dcfSSatish Balay tao->ksp_its+=kspits; 71a7e14dcfSSatish Balay 72a7e14dcfSSatish Balay ierr = VecCopy(tao->stepdirection,ssls->w);CHKERRQ(ierr); 73a7e14dcfSSatish Balay ierr = VecScale(ssls->w,-1.0);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);CHKERRQ(ierr); 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay ierr = VecNorm(ssls->w,NORM_2,&normd);CHKERRQ(ierr); 77a7e14dcfSSatish Balay ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 80a7e14dcfSSatish Balay if (innerd >= -delta*pow(normd, rho)) { 81a7e14dcfSSatish Balay ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr); 82a7e14dcfSSatish Balay ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr); 83a7e14dcfSSatish Balay ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay 86a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 87a7e14dcfSSatish Balay innerd = -innerd; 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0); 90a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr); 91a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 92a7e14dcfSSatish Balay } 93a7e14dcfSSatish Balay PetscFunctionReturn(0); 94a7e14dcfSSatish Balay } 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay #undef __FUNCT__ 97a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_SSFLS" 98441846f8SBarry Smith PetscErrorCode TaoDestroy_SSFLS(Tao tao) 99a7e14dcfSSatish Balay { 100a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 101a7e14dcfSSatish Balay PetscErrorCode ierr; 102a7e14dcfSSatish Balay 103a7e14dcfSSatish Balay PetscFunctionBegin; 104a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr); 105a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->w);CHKERRQ(ierr); 106a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr); 107a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->da);CHKERRQ(ierr); 108a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->db);CHKERRQ(ierr); 109a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr); 110a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr); 111a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 112a7e14dcfSSatish Balay PetscFunctionReturn(0); 113a7e14dcfSSatish Balay } 114a7e14dcfSSatish Balay 115a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 116a7e14dcfSSatish Balay EXTERN_C_BEGIN 117a7e14dcfSSatish Balay #undef __FUNCT__ 118a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_SSFLS" 119441846f8SBarry Smith PetscErrorCode TaoCreate_SSFLS(Tao tao) 120a7e14dcfSSatish Balay { 121a7e14dcfSSatish Balay TAO_SSLS *ssls; 122a7e14dcfSSatish Balay PetscErrorCode ierr; 123a7e14dcfSSatish Balay const char *armijo_type = TAOLINESEARCH_ARMIJO; 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay PetscFunctionBegin; 1263c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr); 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 137a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 138a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr); 139a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch); 14047a47007SBarry Smith /* Linesearch objective and objectivegradient routines are set in solve routine */ 141a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay tao->max_it = 2000; 144a7e14dcfSSatish Balay tao->max_funcs = 4000; 1456f4723b1SBarry Smith tao->fatol = 0; 1466f4723b1SBarry Smith tao->frtol = 0; 1476f4723b1SBarry Smith tao->grtol=0; 1486f4723b1SBarry Smith tao->grtol=0; 1496f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1506f4723b1SBarry Smith tao->gatol = 1.0e-6; 1516f4723b1SBarry Smith tao->fmin = 1.0e-4; 1526f4723b1SBarry Smith #else 153a7e14dcfSSatish Balay tao->gatol = 1.0e-16; 154a7e14dcfSSatish Balay tao->fmin = 1.0e-8; 1556f4723b1SBarry Smith #endif 156a7e14dcfSSatish Balay PetscFunctionReturn(0); 157a7e14dcfSSatish Balay } 158a7e14dcfSSatish Balay EXTERN_C_END 159a7e14dcfSSatish Balay 160