xref: /petsc/src/tao/complementarity/impls/ssls/ssfls.c (revision 3c9e27cfca911a7d7e3219758be42726e83c4ab2)
1a7e14dcfSSatish Balay #include "ssls.h"
2a7e14dcfSSatish Balay 
3a7e14dcfSSatish Balay #undef __FUNCT__
4a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_SSFLS"
5a7e14dcfSSatish Balay PetscErrorCode TaoSetUp_SSFLS(TaoSolver tao)
6a7e14dcfSSatish Balay {
7a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
8a7e14dcfSSatish Balay   PetscErrorCode ierr;
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay   PetscFunctionBegin;
11a7e14dcfSSatish Balay 
12a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->gradient); CHKERRQ(ierr);
13a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->stepdirection); CHKERRQ(ierr);
14a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->w); CHKERRQ(ierr);
15a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->ff); CHKERRQ(ierr);
16a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->dpsi); CHKERRQ(ierr);
17a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->da); CHKERRQ(ierr);
18a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->db); CHKERRQ(ierr);
19a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t1); CHKERRQ(ierr);
20a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t2); CHKERRQ(ierr);
21a7e14dcfSSatish Balay   if (!tao->XL) {
22a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XL); CHKERRQ(ierr);
23a7e14dcfSSatish Balay   }
24a7e14dcfSSatish Balay   if (!tao->XU) {
25a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XU); CHKERRQ(ierr);
26a7e14dcfSSatish Balay   }
27a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU); CHKERRQ(ierr);
28a7e14dcfSSatish Balay   PetscFunctionReturn(0);
29a7e14dcfSSatish Balay }
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay 
32a7e14dcfSSatish Balay #undef __FUNCT__
33a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_SSFLS"
34a7e14dcfSSatish Balay static PetscErrorCode TaoSolve_SSFLS(TaoSolver tao)
35a7e14dcfSSatish Balay {
36a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
37a7e14dcfSSatish Balay   PetscReal psi, ndpsi, normd, innerd, t=0;
38a7e14dcfSSatish Balay   PetscReal delta, rho;
39a7e14dcfSSatish Balay   PetscInt iter=0,kspits;
40a7e14dcfSSatish Balay   TaoSolverTerminationReason reason;
41a7e14dcfSSatish Balay   TaoLineSearchTerminationReason ls_reason;
42a7e14dcfSSatish Balay   PetscErrorCode ierr;
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   PetscFunctionBegin;
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   /* Assume that Setup has been called!
47a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
48a7e14dcfSSatish Balay 
49a7e14dcfSSatish Balay   delta = ssls->delta;
50a7e14dcfSSatish Balay   rho = ssls->rho;
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao); CHKERRQ(ierr);
53a7e14dcfSSatish Balay   /* Project solution inside bounds */
54a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution); CHKERRQ(ierr);
55a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);
56a7e14dcfSSatish Balay   CHKERRQ(ierr);
57a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao); CHKERRQ(ierr);
58a7e14dcfSSatish Balay 
59a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
60a7e14dcfSSatish Balay      current iterate */
61a7e14dcfSSatish Balay   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi); CHKERRQ(ierr);
62a7e14dcfSSatish Balay   ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi); CHKERRQ(ierr);
63a7e14dcfSSatish Balay 
64a7e14dcfSSatish Balay 
65a7e14dcfSSatish Balay   while (1) {
66a7e14dcfSSatish Balay     ierr=PetscInfo3(tao, "iter: %D, merit: %G, ndpsi: %G\n",
67a7e14dcfSSatish Balay 		       iter, ssls->merit, ndpsi); CHKERRQ(ierr);
68a7e14dcfSSatish Balay     /* Check the termination criteria */
69a7e14dcfSSatish Balay     ierr = TaoMonitor(tao,iter++,ssls->merit,ndpsi,0.0,t,&reason);
70a7e14dcfSSatish Balay            CHKERRQ(ierr);
71a7e14dcfSSatish Balay     if (reason!=TAO_CONTINUE_ITERATING) break;
72a7e14dcfSSatish Balay 
73a7e14dcfSSatish Balay     /* Calculate direction.  (Really negative of newton direction.  Therefore,
74a7e14dcfSSatish Balay        rest of the code uses -d.) */
75a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre,ssls->matflag); CHKERRQ(ierr);
76a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection); CHKERRQ(ierr);
77a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&kspits); CHKERRQ(ierr);
78a7e14dcfSSatish Balay     tao->ksp_its+=kspits;
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay     ierr = VecCopy(tao->stepdirection,ssls->w); CHKERRQ(ierr);
81a7e14dcfSSatish Balay     ierr = VecScale(ssls->w,-1.0); CHKERRQ(ierr);
82a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w); CHKERRQ(ierr);
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay     ierr = VecNorm(ssls->w,NORM_2,&normd); CHKERRQ(ierr);
85a7e14dcfSSatish Balay     ierr = VecDot(ssls->w,ssls->dpsi,&innerd); CHKERRQ(ierr);
86a7e14dcfSSatish Balay 
87a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
88a7e14dcfSSatish Balay     if (innerd >= -delta*pow(normd, rho)) {
89a7e14dcfSSatish Balay       ierr = PetscInfo(tao, "newton direction not descent\n"); CHKERRQ(ierr);
90a7e14dcfSSatish Balay       ierr = VecCopy(ssls->dpsi,tao->stepdirection); CHKERRQ(ierr);
91a7e14dcfSSatish Balay       ierr = VecDot(ssls->w,ssls->dpsi,&innerd); CHKERRQ(ierr);
92a7e14dcfSSatish Balay     }
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr);
95a7e14dcfSSatish Balay     innerd = -innerd;
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
98a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason); CHKERRQ(ierr);
99a7e14dcfSSatish Balay     ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi); CHKERRQ(ierr);
100a7e14dcfSSatish Balay   }
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay   PetscFunctionReturn(0);
103a7e14dcfSSatish Balay }
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay #undef __FUNCT__
106a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_SSFLS"
107a7e14dcfSSatish Balay PetscErrorCode TaoDestroy_SSFLS(TaoSolver tao)
108a7e14dcfSSatish Balay {
109a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
110a7e14dcfSSatish Balay   PetscErrorCode ierr;
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay   PetscFunctionBegin;
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->ff); CHKERRQ(ierr);
115a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->w); CHKERRQ(ierr);
116a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dpsi); CHKERRQ(ierr);
117a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->da); CHKERRQ(ierr);
118a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->db); CHKERRQ(ierr);
119a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t1); CHKERRQ(ierr);
120a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t2); CHKERRQ(ierr);
121a7e14dcfSSatish Balay   ierr = PetscFree(tao->data); CHKERRQ(ierr);
122a7e14dcfSSatish Balay   tao->data = PETSC_NULL;
123a7e14dcfSSatish Balay   PetscFunctionReturn(0);
124a7e14dcfSSatish Balay }
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
127a7e14dcfSSatish Balay EXTERN_C_BEGIN
128a7e14dcfSSatish Balay #undef __FUNCT__
129a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_SSFLS"
130a7e14dcfSSatish Balay PetscErrorCode TaoCreate_SSFLS(TaoSolver tao)
131a7e14dcfSSatish Balay {
132a7e14dcfSSatish Balay   TAO_SSLS *ssls;
133a7e14dcfSSatish Balay   PetscErrorCode ierr;
134a7e14dcfSSatish Balay   const char *armijo_type = TAOLINESEARCH_ARMIJO;
135a7e14dcfSSatish Balay 
136a7e14dcfSSatish Balay   PetscFunctionBegin;
137*3c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&ssls); CHKERRQ(ierr);
138a7e14dcfSSatish Balay   tao->data = (void*)ssls;
139a7e14dcfSSatish Balay   tao->ops->solve=TaoSolve_SSFLS;
140a7e14dcfSSatish Balay   tao->ops->setup=TaoSetUp_SSFLS;
141a7e14dcfSSatish Balay   tao->ops->view=TaoView_SSLS;
142a7e14dcfSSatish Balay   tao->ops->setfromoptions=TaoSetFromOptions_SSLS;
143a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_SSFLS;
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   ssls->delta = 1e-10;
146a7e14dcfSSatish Balay   ssls->rho = 2.1;
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch); CHKERRQ(ierr);
149a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,armijo_type); CHKERRQ(ierr);
150a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);
151a7e14dcfSSatish Balay   /* Linesearch objective and objectivegradient routines are
152a7e14dcfSSatish Balay      set in solve routine */
153a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp); CHKERRQ(ierr);
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay   tao->max_it = 2000;
156a7e14dcfSSatish Balay   tao->max_funcs = 4000;
157a7e14dcfSSatish Balay   tao->fatol = 0; tao->frtol = 0; tao->grtol=0; tao->grtol=0;
158a7e14dcfSSatish Balay   tao->gatol = 1.0e-16;
159a7e14dcfSSatish Balay   tao->fmin = 1.0e-8;
160a7e14dcfSSatish Balay 
161a7e14dcfSSatish Balay   PetscFunctionReturn(0);
162a7e14dcfSSatish Balay }
163a7e14dcfSSatish Balay EXTERN_C_END
164a7e14dcfSSatish Balay 
165