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