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