1 #include <../src/tao/complementarity/impls/ssls/ssls.h> 2 3 PetscErrorCode TaoSetUp_SSFLS(Tao tao) 4 { 5 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 6 PetscErrorCode ierr; 7 8 PetscFunctionBegin; 9 ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 10 ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 11 ierr = VecDuplicate(tao->solution,&ssls->w);CHKERRQ(ierr); 12 ierr = VecDuplicate(tao->solution,&ssls->ff);CHKERRQ(ierr); 13 ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr); 14 ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr); 15 ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr); 16 ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr); 17 ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr); 18 if (!tao->XL) { 19 ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr); 20 } 21 if (!tao->XU) { 22 ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr); 23 } 24 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 25 PetscFunctionReturn(0); 26 } 27 28 static PetscErrorCode TaoSolve_SSFLS(Tao tao) 29 { 30 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 31 PetscReal psi, ndpsi, normd, innerd, t=0; 32 PetscReal delta, rho; 33 TaoLineSearchConvergedReason ls_reason; 34 PetscErrorCode ierr; 35 36 PetscFunctionBegin; 37 /* Assume that Setup has been called! 38 Set the structure for the Jacobian and create a linear solver. */ 39 delta = ssls->delta; 40 rho = ssls->rho; 41 42 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 43 /* Project solution inside bounds */ 44 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 45 ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr); 46 ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 47 48 /* Calculate the function value and fischer function value at the 49 current iterate */ 50 ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr); 51 ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 52 53 tao->reason = TAO_CONTINUE_ITERATING; 54 while (PETSC_TRUE) { 55 ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr); 56 /* Check the termination criteria */ 57 ierr = TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its);CHKERRQ(ierr); 58 ierr = TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t);CHKERRQ(ierr); 59 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 60 if (tao->reason!=TAO_CONTINUE_ITERATING) break; 61 tao->niter++; 62 63 /* Calculate direction. (Really negative of newton direction. Therefore, 64 rest of the code uses -d.) */ 65 ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr); 66 ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr); 67 ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr); 68 tao->ksp_tot_its+=tao->ksp_its; 69 70 ierr = VecCopy(tao->stepdirection,ssls->w);CHKERRQ(ierr); 71 ierr = VecScale(ssls->w,-1.0);CHKERRQ(ierr); 72 ierr = VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);CHKERRQ(ierr); 73 74 ierr = VecNorm(ssls->w,NORM_2,&normd);CHKERRQ(ierr); 75 ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 76 77 /* Make sure that we have a descent direction */ 78 if (innerd >= -delta*PetscPowReal(normd, rho)) { 79 ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr); 80 ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr); 81 ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr); 82 } 83 84 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 85 innerd = -innerd; 86 87 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 88 ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr); 89 ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 90 } 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode TaoDestroy_SSFLS(Tao tao) 95 { 96 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 97 PetscErrorCode ierr; 98 99 PetscFunctionBegin; 100 ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr); 101 ierr = VecDestroy(&ssls->w);CHKERRQ(ierr); 102 ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr); 103 ierr = VecDestroy(&ssls->da);CHKERRQ(ierr); 104 ierr = VecDestroy(&ssls->db);CHKERRQ(ierr); 105 ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr); 106 ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr); 107 ierr = PetscFree(tao->data);CHKERRQ(ierr); 108 PetscFunctionReturn(0); 109 } 110 111 /* ---------------------------------------------------------- */ 112 /*MC 113 TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving 114 complementarity constraints 115 116 Options Database Keys: 117 + -tao_ssls_delta - descent test fraction 118 - -tao_ssls_rho - descent test power 119 120 Level: beginner 121 M*/ 122 123 PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao) 124 { 125 TAO_SSLS *ssls; 126 PetscErrorCode ierr; 127 const char *armijo_type = TAOLINESEARCHARMIJO; 128 129 PetscFunctionBegin; 130 ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr); 131 tao->data = (void*)ssls; 132 tao->ops->solve=TaoSolve_SSFLS; 133 tao->ops->setup=TaoSetUp_SSFLS; 134 tao->ops->view=TaoView_SSLS; 135 tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 136 tao->ops->destroy = TaoDestroy_SSFLS; 137 138 ssls->delta = 1e-10; 139 ssls->rho = 2.1; 140 141 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 142 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 143 ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr); 144 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 145 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 146 /* Linesearch objective and objectivegradient routines are set in solve routine */ 147 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 148 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 149 ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 150 151 /* Override default settings (unless already changed) */ 152 if (!tao->max_it_changed) tao->max_it = 2000; 153 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 154 if (!tao->gttol_changed) tao->gttol = 0; 155 if (!tao->grtol_changed) tao->grtol = 0; 156 #if defined(PETSC_USE_REAL_SINGLE) 157 if (!tao->gatol_changed) tao->gatol = 1.0e-6; 158 if (!tao->fmin_changed) tao->fmin = 1.0e-4; 159 #else 160 if (!tao->gatol_changed) tao->gatol = 1.0e-16; 161 if (!tao->fmin_changed) tao->fmin = 1.0e-8; 162 #endif 163 PetscFunctionReturn(0); 164 } 165