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