xref: /petsc/src/tao/complementarity/impls/ssls/ssils.c (revision e4cb33bb7dbdbae9285fba102465ca0f1dcb3977)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
3a7e14dcfSSatish Balay #undef __FUNCT__
4a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_SSILS"
5441846f8SBarry Smith PetscErrorCode TaoSetUp_SSILS(Tao tao)
6a7e14dcfSSatish Balay {
7a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
8a7e14dcfSSatish Balay   PetscErrorCode ierr;
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay   PetscFunctionBegin;
11a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
12a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
13a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->ff);CHKERRQ(ierr);
14a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr);
15a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr);
16a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr);
17a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr);
18a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr);
19a7e14dcfSSatish Balay   PetscFunctionReturn(0);
20a7e14dcfSSatish Balay }
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay #undef __FUNCT__
23a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_SSILS"
24441846f8SBarry Smith PetscErrorCode TaoDestroy_SSILS(Tao tao)
25a7e14dcfSSatish Balay {
26a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
27a7e14dcfSSatish Balay   PetscErrorCode ierr;
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay   PetscFunctionBegin;
30a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr);
31a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr);
32a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->da);CHKERRQ(ierr);
33a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->db);CHKERRQ(ierr);
34a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr);
35a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr);
36a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
37a7e14dcfSSatish Balay   PetscFunctionReturn(0);
38a7e14dcfSSatish Balay }
39a7e14dcfSSatish Balay 
40a7e14dcfSSatish Balay #undef __FUNCT__
41a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_SSILS"
42441846f8SBarry Smith static PetscErrorCode TaoSolve_SSILS(Tao tao)
43a7e14dcfSSatish Balay {
44a7e14dcfSSatish Balay   TAO_SSLS                     *ssls = (TAO_SSLS *)tao->data;
45a7e14dcfSSatish Balay   PetscReal                    psi, ndpsi, normd, innerd, t=0;
46a7e14dcfSSatish Balay   PetscReal                    delta, rho;
47a7e14dcfSSatish Balay   PetscInt                     iter=0,kspits;
48*e4cb33bbSBarry Smith   TaoConvergedReason           reason;
49*e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
50a7e14dcfSSatish Balay   PetscErrorCode               ierr;
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay   PetscFunctionBegin;
53a7e14dcfSSatish Balay   /* Assume that Setup has been called!
54a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
55a7e14dcfSSatish Balay   delta = ssls->delta;
56a7e14dcfSSatish Balay   rho = ssls->rho;
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
59a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
6047a47007SBarry Smith   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr);
61a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr);
62a7e14dcfSSatish Balay 
63a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
64a7e14dcfSSatish Balay      current iterate */
65a7e14dcfSSatish Balay   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr);
66a7e14dcfSSatish Balay   ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay   while (1) {
6947a47007SBarry Smith     ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",iter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr);
70a7e14dcfSSatish Balay     /* Check the termination criteria */
7147a47007SBarry Smith     ierr = TaoMonitor(tao,iter++,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr);
72a7e14dcfSSatish Balay     if (reason!=TAO_CONTINUE_ITERATING) break;
73a7e14dcfSSatish Balay 
74a7e14dcfSSatish Balay     /* Calculate direction.  (Really negative of newton direction.  Therefore,
75a7e14dcfSSatish Balay        rest of the code uses -d.) */
76a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre,ssls->matflag);CHKERRQ(ierr);
77a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr);
78a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
79a7e14dcfSSatish Balay     tao->ksp_its+=kspits;
80a7e14dcfSSatish Balay     ierr = VecNorm(tao->stepdirection,NORM_2,&normd);CHKERRQ(ierr);
81a7e14dcfSSatish Balay     ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr);
82a7e14dcfSSatish Balay 
83a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
84a7e14dcfSSatish Balay     if (innerd <= delta*pow(normd, rho)) {
85a7e14dcfSSatish Balay       ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr);
86a7e14dcfSSatish Balay       ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr);
87a7e14dcfSSatish Balay       ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr);
88a7e14dcfSSatish Balay     }
89a7e14dcfSSatish Balay 
90a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
91a7e14dcfSSatish Balay     innerd = -innerd;
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
94a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr);
95a7e14dcfSSatish Balay     ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
96a7e14dcfSSatish Balay   }
97a7e14dcfSSatish Balay   PetscFunctionReturn(0);
98a7e14dcfSSatish Balay }
99a7e14dcfSSatish Balay 
100a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
101a7e14dcfSSatish Balay EXTERN_C_BEGIN
102a7e14dcfSSatish Balay #undef __FUNCT__
103a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_SSILS"
104441846f8SBarry Smith PetscErrorCode TaoCreate_SSILS(Tao tao)
105a7e14dcfSSatish Balay {
106a7e14dcfSSatish Balay   TAO_SSLS       *ssls;
107a7e14dcfSSatish Balay   PetscErrorCode ierr;
108a7e14dcfSSatish Balay   const char     *armijo_type = TAOLINESEARCH_ARMIJO;
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay   PetscFunctionBegin;
1113c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr);
112a7e14dcfSSatish Balay   tao->data = (void*)ssls;
113a7e14dcfSSatish Balay   tao->ops->solve=TaoSolve_SSILS;
114a7e14dcfSSatish Balay   tao->ops->setup=TaoSetUp_SSILS;
115a7e14dcfSSatish Balay   tao->ops->view=TaoView_SSLS;
116a7e14dcfSSatish Balay   tao->ops->setfromoptions=TaoSetFromOptions_SSLS;
117a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_SSILS;
118a7e14dcfSSatish Balay 
119a7e14dcfSSatish Balay   ssls->delta = 1e-10;
120a7e14dcfSSatish Balay   ssls->rho = 2.1;
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
123a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr);
124a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);
125a7e14dcfSSatish Balay   /* Note: linesearch objective and objectivegradient routines are set in solve routine */
126a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay   tao->max_it = 2000;
129a7e14dcfSSatish Balay   tao->max_funcs = 4000;
1306f4723b1SBarry Smith   tao->fatol = 0;
1316f4723b1SBarry Smith   tao->frtol = 0;
1326f4723b1SBarry Smith   tao->gttol=0;
1336f4723b1SBarry Smith   tao->grtol=0;
1346f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1356f4723b1SBarry Smith   tao->gatol = 1.0e-6;
1366f4723b1SBarry Smith   tao->fmin = 1.0e-4;
1376f4723b1SBarry Smith #else
138a7e14dcfSSatish Balay   tao->gatol = 1.0e-16;
139a7e14dcfSSatish Balay   tao->fmin = 1.0e-8;
1406f4723b1SBarry Smith #endif
141a7e14dcfSSatish Balay   PetscFunctionReturn(0);
142a7e14dcfSSatish Balay }
143a7e14dcfSSatish Balay EXTERN_C_END
144a7e14dcfSSatish Balay 
145