Actual source code: ssils.c
petsc-3.8.4 2018-03-24
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: TaoConvergedReason reason;
42: TaoLineSearchConvergedReason ls_reason;
43: PetscErrorCode ierr;
46: /* Assume that Setup has been called!
47: Set the structure for the Jacobian and create a linear solver. */
48: delta = ssls->delta;
49: rho = ssls->rho;
51: TaoComputeVariableBounds(tao);
52: VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
53: TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);
54: TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);
56: /* Calculate the function value and fischer function value at the
57: current iterate */
58: TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);
59: VecNorm(ssls->dpsi,NORM_2,&ndpsi);
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: TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t,&reason);
65: if (reason!=TAO_CONTINUE_ITERATING) break;
66: tao->niter++;
68: /* Calculate direction. (Really negative of newton direction. Therefore,
69: rest of the code uses -d.) */
70: KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);
71: KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);
72: KSPGetIterationNumber(tao->ksp,&tao->ksp_its);
73: tao->ksp_tot_its+=tao->ksp_its;
74: VecNorm(tao->stepdirection,NORM_2,&normd);
75: VecDot(tao->stepdirection,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(tao->stepdirection,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: /* ---------------------------------------------------------- */
95: /*MC
96: TAOSSILS - semi-smooth infeasible linesearch algorithm for solving
97: complementarity constraints
99: Options Database Keys:
100: + -tao_ssls_delta - descent test fraction
101: - -tao_ssls_rho - descent test power
103: Level: beginner
104: M*/
105: PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao)
106: {
107: TAO_SSLS *ssls;
109: const char *armijo_type = TAOLINESEARCHARMIJO;
112: PetscNewLog(tao,&ssls);
113: tao->data = (void*)ssls;
114: tao->ops->solve=TaoSolve_SSILS;
115: tao->ops->setup=TaoSetUp_SSILS;
116: tao->ops->view=TaoView_SSLS;
117: tao->ops->setfromoptions=TaoSetFromOptions_SSLS;
118: tao->ops->destroy = TaoDestroy_SSILS;
120: ssls->delta = 1e-10;
121: ssls->rho = 2.1;
123: TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);
124: TaoLineSearchSetType(tao->linesearch,armijo_type);
125: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
126: TaoLineSearchSetFromOptions(tao->linesearch);
127: /* Note: linesearch objective and objectivegradient routines are set in solve routine */
128: KSPCreate(((PetscObject)tao)->comm,&tao->ksp);
129: KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);
131: /* Override default settings (unless already changed) */
132: if (!tao->max_it_changed) tao->max_it = 2000;
133: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
134: if (!tao->gttol_changed) tao->gttol = 0;
135: if (!tao->grtol_changed) tao->grtol = 0;
136: #if defined(PETSC_USE_REAL_SINGLE)
137: if (!tao->gatol_changed) tao->gatol = 1.0e-6;
138: if (!tao->fmin_changed) tao->fmin = 1.0e-4;
139: #else
140: if (!tao->gatol_changed) tao->gatol = 1.0e-16;
141: if (!tao->fmin_changed) tao->fmin = 1.0e-8;
142: #endif
143: return(0);
144: }