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