xref: /petsc/src/tao/complementarity/impls/ssls/ssils.c (revision 63b1541575edd3f8a34e081896782a8d7a30de53)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
3441846f8SBarry Smith PetscErrorCode TaoSetUp_SSILS(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->ff);CHKERRQ(ierr);
12a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr);
13a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr);
14a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr);
15a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr);
16a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr);
17a7e14dcfSSatish Balay   PetscFunctionReturn(0);
18a7e14dcfSSatish Balay }
19a7e14dcfSSatish Balay 
20441846f8SBarry Smith PetscErrorCode TaoDestroy_SSILS(Tao tao)
21a7e14dcfSSatish Balay {
22a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
23a7e14dcfSSatish Balay   PetscErrorCode ierr;
24a7e14dcfSSatish Balay 
25a7e14dcfSSatish Balay   PetscFunctionBegin;
26a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr);
27a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr);
28a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->da);CHKERRQ(ierr);
29a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->db);CHKERRQ(ierr);
30a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr);
31a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr);
32a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
33a7e14dcfSSatish Balay   PetscFunctionReturn(0);
34a7e14dcfSSatish Balay }
35a7e14dcfSSatish Balay 
36441846f8SBarry Smith static PetscErrorCode TaoSolve_SSILS(Tao tao)
37a7e14dcfSSatish Balay {
38a7e14dcfSSatish Balay   TAO_SSLS                     *ssls = (TAO_SSLS *)tao->data;
39a7e14dcfSSatish Balay   PetscReal                    psi, ndpsi, normd, innerd, t=0;
40a7e14dcfSSatish Balay   PetscReal                    delta, rho;
41e4cb33bbSBarry Smith   TaoConvergedReason           reason;
42e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
43a7e14dcfSSatish Balay   PetscErrorCode               ierr;
44a7e14dcfSSatish Balay 
45a7e14dcfSSatish Balay   PetscFunctionBegin;
46a7e14dcfSSatish Balay   /* Assume that Setup has been called!
47a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
48a7e14dcfSSatish Balay   delta = ssls->delta;
49a7e14dcfSSatish Balay   rho = ssls->rho;
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
52a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
5347a47007SBarry Smith   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr);
54a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr);
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
57a7e14dcfSSatish Balay      current iterate */
58a7e14dcfSSatish Balay   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr);
59a7e14dcfSSatish Balay   ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
60a7e14dcfSSatish Balay 
61d90ca52eSBarry Smith   while (PETSC_TRUE) {
628931d482SJason Sarich     ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr);
63a7e14dcfSSatish Balay     /* Check the termination criteria */
64e6d4cb7fSJason Sarich     ierr = TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr);
65a7e14dcfSSatish Balay     if (reason!=TAO_CONTINUE_ITERATING) break;
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     ierr = VecNorm(tao->stepdirection,NORM_2,&normd);CHKERRQ(ierr);
75a7e14dcfSSatish Balay     ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr);
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
781118d4bcSLisandro Dalcin     if (innerd <= delta*PetscPowReal(normd, rho)) {
79a7e14dcfSSatish Balay       ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr);
80a7e14dcfSSatish Balay       ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr);
81a7e14dcfSSatish Balay       ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr);
82a7e14dcfSSatish Balay     }
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
85a7e14dcfSSatish Balay     innerd = -innerd;
86a7e14dcfSSatish Balay 
87302440fdSBarry Smith     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
88a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr);
89a7e14dcfSSatish Balay     ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
90a7e14dcfSSatish Balay   }
91a7e14dcfSSatish Balay   PetscFunctionReturn(0);
92a7e14dcfSSatish Balay }
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
951522df2eSJason Sarich /*MC
961522df2eSJason Sarich    TAOSSILS - semi-smooth infeasible linesearch algorithm for solving
971522df2eSJason Sarich        complementarity constraints
981522df2eSJason Sarich 
991522df2eSJason Sarich    Options Database Keys:
1001522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
1011522df2eSJason Sarich - -tao_ssls_rho - descent test power
1021522df2eSJason Sarich 
1031eb8069cSJason Sarich    Level: beginner
1041522df2eSJason Sarich M*/
105728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao)
106a7e14dcfSSatish Balay {
107a7e14dcfSSatish Balay   TAO_SSLS       *ssls;
108a7e14dcfSSatish Balay   PetscErrorCode ierr;
1098caf6e8cSBarry Smith   const char     *armijo_type = TAOLINESEARCHARMIJO;
110a7e14dcfSSatish Balay 
111a7e14dcfSSatish Balay   PetscFunctionBegin;
1123c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr);
113a7e14dcfSSatish Balay   tao->data = (void*)ssls;
114a7e14dcfSSatish Balay   tao->ops->solve=TaoSolve_SSILS;
115a7e14dcfSSatish Balay   tao->ops->setup=TaoSetUp_SSILS;
116a7e14dcfSSatish Balay   tao->ops->view=TaoView_SSLS;
117a7e14dcfSSatish Balay   tao->ops->setfromoptions=TaoSetFromOptions_SSLS;
118a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_SSILS;
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay   ssls->delta = 1e-10;
121a7e14dcfSSatish Balay   ssls->rho = 2.1;
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
124*63b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
125a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr);
1265d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
127302440fdSBarry Smith   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
128a7e14dcfSSatish Balay   /* Note: linesearch objective and objectivegradient routines are set in solve routine */
129a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
130*63b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1315d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
132a7e14dcfSSatish Balay 
1336552cf8aSJason Sarich   /* Override default settings (unless already changed) */
1346552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
1356552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
1366552cf8aSJason Sarich   if (!tao->gttol_changed) tao->gttol = 0;
1376552cf8aSJason Sarich   if (!tao->grtol_changed) tao->grtol = 0;
1386f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1396552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
1406552cf8aSJason Sarich   if (!tao->fmin_changed)  tao->fmin = 1.0e-4;
1416f4723b1SBarry Smith #else
1426552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
1436552cf8aSJason Sarich   if (!tao->fmin_changed)  tao->fmin = 1.0e-8;
1446f4723b1SBarry Smith #endif
145a7e14dcfSSatish Balay   PetscFunctionReturn(0);
146a7e14dcfSSatish Balay }
147