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