Actual source code: ssfls.c
petsc-3.6.4 2016-04-12
1: #include <../src/tao/complementarity/impls/ssls/ssls.h>
5: PetscErrorCode TaoSetUp_SSFLS(Tao tao)
6: {
7: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
11: VecDuplicate(tao->solution,&tao->gradient);
12: VecDuplicate(tao->solution,&tao->stepdirection);
13: VecDuplicate(tao->solution,&ssls->w);
14: VecDuplicate(tao->solution,&ssls->ff);
15: VecDuplicate(tao->solution,&ssls->dpsi);
16: VecDuplicate(tao->solution,&ssls->da);
17: VecDuplicate(tao->solution,&ssls->db);
18: VecDuplicate(tao->solution,&ssls->t1);
19: VecDuplicate(tao->solution,&ssls->t2);
20: if (!tao->XL) {
21: VecDuplicate(tao->solution,&tao->XL);
22: }
23: if (!tao->XU) {
24: VecDuplicate(tao->solution,&tao->XU);
25: }
26: TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
27: return(0);
28: }
32: static PetscErrorCode TaoSolve_SSFLS(Tao tao)
33: {
34: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
35: PetscReal psi, ndpsi, normd, innerd, t=0;
36: PetscReal delta, rho;
37: TaoConvergedReason reason;
38: TaoLineSearchConvergedReason ls_reason;
39: PetscErrorCode ierr;
42: /* Assume that Setup has been called!
43: Set the structure for the Jacobian and create a linear solver. */
44: delta = ssls->delta;
45: rho = ssls->rho;
47: TaoComputeVariableBounds(tao);
48: /* Project solution inside bounds */
49: VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
50: TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);
51: TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);
53: /* Calculate the function value and fischer function value at the
54: current iterate */
55: TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);
56: VecNorm(ssls->dpsi,NORM_2,&ndpsi);
58: while (PETSC_TRUE) {
59: ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);
60: /* Check the termination criteria */
61: TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t,&reason);
62: if (reason!=TAO_CONTINUE_ITERATING) break;
63: tao->niter++;
65: /* Calculate direction. (Really negative of newton direction. Therefore,
66: rest of the code uses -d.) */
67: KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);
68: KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);
69: KSPGetIterationNumber(tao->ksp,&tao->ksp_its);
70: tao->ksp_tot_its+=tao->ksp_its;
72: VecCopy(tao->stepdirection,ssls->w);
73: VecScale(ssls->w,-1.0);
74: VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);
76: VecNorm(ssls->w,NORM_2,&normd);
77: VecDot(ssls->w,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(ssls->w,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: }
98: PetscErrorCode TaoDestroy_SSFLS(Tao tao)
99: {
100: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
104: VecDestroy(&ssls->ff);
105: VecDestroy(&ssls->w);
106: VecDestroy(&ssls->dpsi);
107: VecDestroy(&ssls->da);
108: VecDestroy(&ssls->db);
109: VecDestroy(&ssls->t1);
110: VecDestroy(&ssls->t2);
111: PetscFree(tao->data);
112: return(0);
113: }
115: /* ---------------------------------------------------------- */
116: /*MC
117: TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving
118: complementarity constraints
120: Options Database Keys:
121: + -tao_ssls_delta - descent test fraction
122: - -tao_ssls_rho - descent test power
124: Level: beginner
125: M*/
129: PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao)
130: {
131: TAO_SSLS *ssls;
133: const char *armijo_type = TAOLINESEARCHARMIJO;
136: PetscNewLog(tao,&ssls);
137: tao->data = (void*)ssls;
138: tao->ops->solve=TaoSolve_SSFLS;
139: tao->ops->setup=TaoSetUp_SSFLS;
140: tao->ops->view=TaoView_SSLS;
141: tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
142: tao->ops->destroy = TaoDestroy_SSFLS;
144: ssls->delta = 1e-10;
145: ssls->rho = 2.1;
147: TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);
148: TaoLineSearchSetType(tao->linesearch,armijo_type);
149: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
150: TaoLineSearchSetFromOptions(tao->linesearch);
151: /* Linesearch objective and objectivegradient routines are set in solve routine */
152: KSPCreate(((PetscObject)tao)->comm,&tao->ksp);
153: KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);
155: /* Override default settings (unless already changed) */
156: if (!tao->max_it_changed) tao->max_it = 2000;
157: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
158: if (!tao->fatol_changed) tao->fatol = 0;
159: if (!tao->frtol_changed) tao->frtol = 0;
160: if (!tao->gttol_changed) tao->gttol = 0;
161: if (!tao->grtol_changed) tao->grtol = 0;
162: #if defined(PETSC_USE_REAL_SINGLE)
163: if (!tao->gatol_changed) tao->gatol = 1.0e-6;
164: if (!tao->fmin_changed) tao->fmin = 1.0e-4;
165: #else
166: if (!tao->gatol_changed) tao->gatol = 1.0e-16;
167: if (!tao->fmin_changed) tao->fmin = 1.0e-8;
168: #endif
169: return(0);
170: }