xref: /petsc/src/tao/complementarity/impls/ssls/ssfls.c (revision e1e80dc898c3d5fa028e909e947004e745fb92d9)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
3441846f8SBarry Smith PetscErrorCode TaoSetUp_SSFLS(Tao tao)
4a7e14dcfSSatish Balay {
5a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
6a7e14dcfSSatish Balay   PetscErrorCode ierr;
7a7e14dcfSSatish Balay 
8a7e14dcfSSatish Balay   PetscFunctionBegin;
9a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
10a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
11a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->w);CHKERRQ(ierr);
12a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->ff);CHKERRQ(ierr);
13a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr);
14a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr);
15a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr);
16a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr);
17a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr);
18a7e14dcfSSatish Balay   if (!tao->XL) {
19a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
20a7e14dcfSSatish Balay   }
21a7e14dcfSSatish Balay   if (!tao->XU) {
22a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
23a7e14dcfSSatish Balay   }
24a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
25a7e14dcfSSatish Balay   PetscFunctionReturn(0);
26a7e14dcfSSatish Balay }
27a7e14dcfSSatish Balay 
28441846f8SBarry Smith static PetscErrorCode TaoSolve_SSFLS(Tao tao)
29a7e14dcfSSatish Balay {
30a7e14dcfSSatish Balay   TAO_SSLS                     *ssls = (TAO_SSLS *)tao->data;
31a7e14dcfSSatish Balay   PetscReal                    psi, ndpsi, normd, innerd, t=0;
32a7e14dcfSSatish Balay   PetscReal                    delta, rho;
33e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
34a7e14dcfSSatish Balay   PetscErrorCode               ierr;
35a7e14dcfSSatish Balay 
36a7e14dcfSSatish Balay   PetscFunctionBegin;
37a7e14dcfSSatish Balay   /* Assume that Setup has been called!
38a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
39a7e14dcfSSatish Balay   delta = ssls->delta;
40a7e14dcfSSatish Balay   rho = ssls->rho;
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
43a7e14dcfSSatish Balay   /* Project solution inside bounds */
44a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
4547a47007SBarry Smith   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr);
46a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr);
47a7e14dcfSSatish Balay 
48a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
49a7e14dcfSSatish Balay      current iterate */
50a7e14dcfSSatish Balay   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr);
51a7e14dcfSSatish Balay   ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
52a7e14dcfSSatish Balay 
53763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
54d90ca52eSBarry Smith   while (PETSC_TRUE) {
558931d482SJason Sarich     ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr);
56a7e14dcfSSatish Balay     /* Check the termination criteria */
57763847b4SAlp Dener     ierr = TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its);CHKERRQ(ierr);
58763847b4SAlp Dener     ierr = TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t);CHKERRQ(ierr);
59763847b4SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
60763847b4SAlp Dener     if (tao->reason!=TAO_CONTINUE_ITERATING) break;
61*e1e80dc8SAlp Dener 
62*e1e80dc8SAlp Dener     /* Call general purpose update function */
63*e1e80dc8SAlp Dener     if (tao->ops->update) {
64*e1e80dc8SAlp Dener       ierr = (*tao->ops->update)(tao, tao->niter);CHKERRQ(ierr);
65*e1e80dc8SAlp Dener     }
66e6d4cb7fSJason Sarich     tao->niter++;
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay     /* Calculate direction.  (Really negative of newton direction.  Therefore,
69a7e14dcfSSatish Balay        rest of the code uses -d.) */
7023ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr);
71a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr);
72b0026674SJason Sarich     ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr);
73b0026674SJason Sarich     tao->ksp_tot_its+=tao->ksp_its;
74a7e14dcfSSatish Balay 
75a7e14dcfSSatish Balay     ierr = VecCopy(tao->stepdirection,ssls->w);CHKERRQ(ierr);
76a7e14dcfSSatish Balay     ierr = VecScale(ssls->w,-1.0);CHKERRQ(ierr);
77a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w);CHKERRQ(ierr);
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay     ierr = VecNorm(ssls->w,NORM_2,&normd);CHKERRQ(ierr);
80a7e14dcfSSatish Balay     ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr);
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
83d90ca52eSBarry Smith     if (innerd >= -delta*PetscPowReal(normd, rho)) {
84a7e14dcfSSatish Balay       ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr);
85a7e14dcfSSatish Balay       ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr);
86a7e14dcfSSatish Balay       ierr = VecDot(ssls->w,ssls->dpsi,&innerd);CHKERRQ(ierr);
87a7e14dcfSSatish Balay     }
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
90a7e14dcfSSatish Balay     innerd = -innerd;
91a7e14dcfSSatish Balay 
92302440fdSBarry Smith     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
93a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr);
94a7e14dcfSSatish Balay     ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
95a7e14dcfSSatish Balay   }
96a7e14dcfSSatish Balay   PetscFunctionReturn(0);
97a7e14dcfSSatish Balay }
98a7e14dcfSSatish Balay 
99441846f8SBarry Smith PetscErrorCode TaoDestroy_SSFLS(Tao tao)
100a7e14dcfSSatish Balay {
101a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
102a7e14dcfSSatish Balay   PetscErrorCode ierr;
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay   PetscFunctionBegin;
105a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr);
106a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->w);CHKERRQ(ierr);
107a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr);
108a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->da);CHKERRQ(ierr);
109a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->db);CHKERRQ(ierr);
110a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr);
111a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr);
112a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
113a7e14dcfSSatish Balay   PetscFunctionReturn(0);
114a7e14dcfSSatish Balay }
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
1171522df2eSJason Sarich /*MC
1181522df2eSJason Sarich    TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving
1191522df2eSJason Sarich        complementarity constraints
1201522df2eSJason Sarich 
1211522df2eSJason Sarich    Options Database Keys:
1221522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
1231522df2eSJason Sarich - -tao_ssls_rho - descent test power
1241522df2eSJason Sarich 
1251eb8069cSJason Sarich    Level: beginner
1261522df2eSJason Sarich M*/
1271eb8069cSJason Sarich 
128728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao)
129a7e14dcfSSatish Balay {
130a7e14dcfSSatish Balay   TAO_SSLS       *ssls;
131a7e14dcfSSatish Balay   PetscErrorCode ierr;
1328caf6e8cSBarry Smith   const char     *armijo_type = TAOLINESEARCHARMIJO;
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay   PetscFunctionBegin;
1353c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr);
136a7e14dcfSSatish Balay   tao->data = (void*)ssls;
137a7e14dcfSSatish Balay   tao->ops->solve=TaoSolve_SSFLS;
138a7e14dcfSSatish Balay   tao->ops->setup=TaoSetUp_SSFLS;
139a7e14dcfSSatish Balay   tao->ops->view=TaoView_SSLS;
140a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
141a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_SSFLS;
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay   ssls->delta = 1e-10;
144a7e14dcfSSatish Balay   ssls->rho = 2.1;
145a7e14dcfSSatish Balay 
146a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
14763b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
148a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr);
1495d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
150302440fdSBarry Smith   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
15147a47007SBarry Smith   /* Linesearch objective and objectivegradient routines are  set in solve routine */
152a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
15363b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1545d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
155a7e14dcfSSatish Balay 
1566552cf8aSJason Sarich   /* Override default settings (unless already changed) */
1576552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
1586552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
1596552cf8aSJason Sarich   if (!tao->gttol_changed) tao->gttol = 0;
1606552cf8aSJason Sarich   if (!tao->grtol_changed) tao->grtol = 0;
1616f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1626552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
1636552cf8aSJason Sarich   if (!tao->fmin_changed)  tao->fmin = 1.0e-4;
1646f4723b1SBarry Smith #else
1656552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
1666552cf8aSJason Sarich   if (!tao->fmin_changed)  tao->fmin = 1.0e-8;
1676f4723b1SBarry Smith #endif
168a7e14dcfSSatish Balay   PetscFunctionReturn(0);
169a7e14dcfSSatish Balay }
170