Actual source code: ssfls.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <../src/tao/complementarity/impls/ssls/ssls.h>

  3: PetscErrorCode TaoSetUp_SSFLS(Tao tao)
  4: {
  5:   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;

  9:   VecDuplicate(tao->solution,&tao->gradient);
 10:   VecDuplicate(tao->solution,&tao->stepdirection);
 11:   VecDuplicate(tao->solution,&ssls->w);
 12:   VecDuplicate(tao->solution,&ssls->ff);
 13:   VecDuplicate(tao->solution,&ssls->dpsi);
 14:   VecDuplicate(tao->solution,&ssls->da);
 15:   VecDuplicate(tao->solution,&ssls->db);
 16:   VecDuplicate(tao->solution,&ssls->t1);
 17:   VecDuplicate(tao->solution,&ssls->t2);
 18:   if (!tao->XL) {
 19:     VecDuplicate(tao->solution,&tao->XL);
 20:   }
 21:   if (!tao->XU) {
 22:     VecDuplicate(tao->solution,&tao->XU);
 23:   }
 24:   TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
 25:   return(0);
 26: }

 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;

 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;

 42:   TaoComputeVariableBounds(tao);
 43:   /* Project solution inside bounds */
 44:   VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
 45:   TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);
 46:   TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);

 48:   /* Calculate the function value and fischer function value at the
 49:      current iterate */
 50:   TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);
 51:   VecNorm(ssls->dpsi,NORM_2,&ndpsi);

 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);
 56:     /* Check the termination criteria */
 57:     TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its);
 58:     TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t);
 59:     (*tao->ops->convergencetest)(tao,tao->cnvP);
 60:     if (tao->reason!=TAO_CONTINUE_ITERATING) break;
 61:     tao->niter++;

 63:     /* Calculate direction.  (Really negative of newton direction.  Therefore,
 64:        rest of the code uses -d.) */
 65:     KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);
 66:     KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);
 67:     KSPGetIterationNumber(tao->ksp,&tao->ksp_its);
 68:     tao->ksp_tot_its+=tao->ksp_its;

 70:     VecCopy(tao->stepdirection,ssls->w);
 71:     VecScale(ssls->w,-1.0);
 72:     VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);

 74:     VecNorm(ssls->w,NORM_2,&normd);
 75:     VecDot(ssls->w,ssls->dpsi,&innerd);

 77:     /* Make sure that we have a descent direction */
 78:     if (innerd >= -delta*PetscPowReal(normd, rho)) {
 79:       PetscInfo(tao, "newton direction not descent\n");
 80:       VecCopy(ssls->dpsi,tao->stepdirection);
 81:       VecDot(ssls->w,ssls->dpsi,&innerd);
 82:     }

 84:     VecScale(tao->stepdirection, -1.0);
 85:     innerd = -innerd;

 87:     TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
 88:     TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);
 89:     VecNorm(ssls->dpsi,NORM_2,&ndpsi);
 90:   }
 91:   return(0);
 92: }

 94: PetscErrorCode TaoDestroy_SSFLS(Tao tao)
 95: {
 96:   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;

100:   VecDestroy(&ssls->ff);
101:   VecDestroy(&ssls->w);
102:   VecDestroy(&ssls->dpsi);
103:   VecDestroy(&ssls->da);
104:   VecDestroy(&ssls->db);
105:   VecDestroy(&ssls->t1);
106:   VecDestroy(&ssls->t2);
107:   PetscFree(tao->data);
108:   return(0);
109: }

111: /* ---------------------------------------------------------- */
112: /*MC
113:    TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving
114:        complementarity constraints

116:    Options Database Keys:
117: + -tao_ssls_delta - descent test fraction
118: - -tao_ssls_rho - descent test power

120:    Level: beginner
121: M*/

123: PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao)
124: {
125:   TAO_SSLS       *ssls;
127:   const char     *armijo_type = TAOLINESEARCHARMIJO;

130:   PetscNewLog(tao,&ssls);
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;

138:   ssls->delta = 1e-10;
139:   ssls->rho = 2.1;

141:   TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);
142:   PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
143:   TaoLineSearchSetType(tao->linesearch,armijo_type);
144:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
145:   TaoLineSearchSetFromOptions(tao->linesearch);
146:   /* Linesearch objective and objectivegradient routines are  set in solve routine */
147:   KSPCreate(((PetscObject)tao)->comm,&tao->ksp);
148:   PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);
149:   KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);

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:   return(0);
164: }