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