xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision ffad99011bdf8bdff5e8540ef3c49b4fd8d6e6bb)
1aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/ntr/ntr.h>
3a7e14dcfSSatish Balay 
4aaa7dc30SBarry Smith #include <petscksp.h>
5aaa7dc30SBarry Smith #include <petscpc.h>
6aaa7dc30SBarry Smith #include <petsc-private/kspimpl.h>
7aaa7dc30SBarry Smith #include <petsc-private/pcimpl.h>
8a7e14dcfSSatish Balay 
9a7e14dcfSSatish Balay #define NTR_KSP_NASH    0
10a7e14dcfSSatish Balay #define NTR_KSP_STCG    1
11a7e14dcfSSatish Balay #define NTR_KSP_GLTR    2
12a7e14dcfSSatish Balay #define NTR_KSP_TYPES   3
13a7e14dcfSSatish Balay 
14a7e14dcfSSatish Balay #define NTR_PC_NONE     0
15a7e14dcfSSatish Balay #define NTR_PC_AHESS    1
16a7e14dcfSSatish Balay #define NTR_PC_BFGS     2
17a7e14dcfSSatish Balay #define NTR_PC_PETSC    3
18a7e14dcfSSatish Balay #define NTR_PC_TYPES    4
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS   0
21a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS    1
22a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES   2
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay #define NTR_INIT_CONSTANT         0
25a7e14dcfSSatish Balay #define NTR_INIT_DIRECTION        1
26a7e14dcfSSatish Balay #define NTR_INIT_INTERPOLATION    2
27a7e14dcfSSatish Balay #define NTR_INIT_TYPES            3
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay #define NTR_UPDATE_REDUCTION      0
30a7e14dcfSSatish Balay #define NTR_UPDATE_INTERPOLATION  1
31a7e14dcfSSatish Balay #define NTR_UPDATE_TYPES          2
32a7e14dcfSSatish Balay 
3353506e15SBarry Smith static const char *NTR_KSP[64] = {  "nash", "stcg", "gltr"};
34a7e14dcfSSatish Balay 
3553506e15SBarry Smith static const char *NTR_PC[64] = {  "none", "ahess", "bfgs", "petsc"};
36a7e14dcfSSatish Balay 
3753506e15SBarry Smith static const char *BFGS_SCALE[64] = {  "ahess", "bfgs"};
38a7e14dcfSSatish Balay 
3953506e15SBarry Smith static const char *NTR_INIT[64] = {  "constant", "direction", "interpolation"};
40a7e14dcfSSatish Balay 
4153506e15SBarry Smith static const char *NTR_UPDATE[64] = {  "reduction", "interpolation"};
42a7e14dcfSSatish Balay 
43a7e14dcfSSatish Balay /*  Routine for BFGS preconditioner */
44a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec xin, Vec xout);
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay /*
47a7e14dcfSSatish Balay    TaoSolve_NTR - Implements Newton's Method with a trust region approach
48a7e14dcfSSatish Balay    for solving unconstrained minimization problems.
49a7e14dcfSSatish Balay 
50a7e14dcfSSatish Balay    The basic algorithm is taken from MINPACK-2 (dstrn).
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay    TaoSolve_NTR computes a local minimizer of a twice differentiable function
53a7e14dcfSSatish Balay    f by applying a trust region variant of Newton's method.  At each stage
54a7e14dcfSSatish Balay    of the algorithm, we use the prconditioned conjugate gradient method to
55a7e14dcfSSatish Balay    determine an approximate minimizer of the quadratic equation
56a7e14dcfSSatish Balay 
57a7e14dcfSSatish Balay         q(s) = <s, Hs + g>
58a7e14dcfSSatish Balay 
59a7e14dcfSSatish Balay    subject to the trust region constraint
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay         || s ||_M <= radius,
62a7e14dcfSSatish Balay 
63a7e14dcfSSatish Balay    where radius is the trust region radius and M is a symmetric positive
64a7e14dcfSSatish Balay    definite matrix (the preconditioner).  Here g is the gradient and H
65a7e14dcfSSatish Balay    is the Hessian matrix.
66a7e14dcfSSatish Balay 
67a7e14dcfSSatish Balay    Note:  TaoSolve_NTR MUST use the iterative solver KSPNASH, KSPSTCG,
68a7e14dcfSSatish Balay           or KSPGLTR.  Thus, we set KSPNASH, KSPSTCG, or KSPGLTR in this
69a7e14dcfSSatish Balay           routine regardless of what the user may have previously specified.
70a7e14dcfSSatish Balay */
71a7e14dcfSSatish Balay #undef __FUNCT__
72a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NTR"
73441846f8SBarry Smith static PetscErrorCode TaoSolve_NTR(Tao tao)
74a7e14dcfSSatish Balay {
75a7e14dcfSSatish Balay   TAO_NTR            *tr = (TAO_NTR *)tao->data;
76a7e14dcfSSatish Balay   PC                 pc;
77a7e14dcfSSatish Balay   KSPConvergedReason ksp_reason;
78e4cb33bbSBarry Smith   TaoConvergedReason reason;
79a7e14dcfSSatish Balay   PetscReal          fmin, ftrial, prered, actred, kappa, sigma, beta;
80a7e14dcfSSatish Balay   PetscReal          tau, tau_1, tau_2, tau_max, tau_min, max_radius;
81a7e14dcfSSatish Balay   PetscReal          f, gnorm;
82a7e14dcfSSatish Balay 
83a7e14dcfSSatish Balay   PetscReal          delta;
84a7e14dcfSSatish Balay   PetscReal          norm_d;
85a7e14dcfSSatish Balay   PetscErrorCode     ierr;
86a7e14dcfSSatish Balay   PetscInt           iter = 0;
87a7e14dcfSSatish Balay   PetscInt           bfgsUpdates = 0;
88a7e14dcfSSatish Balay   PetscInt           needH;
89a7e14dcfSSatish Balay 
90a7e14dcfSSatish Balay   PetscInt           i_max = 5;
91a7e14dcfSSatish Balay   PetscInt           j_max = 1;
92a7e14dcfSSatish Balay   PetscInt           i, j, N, n, its;
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   PetscFunctionBegin;
95a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
96a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr);
97a7e14dcfSSatish Balay   }
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay   tao->trust = tao->trust0;
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay   /* Modify the radius if it is too large or small */
102a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tr->min_radius);
103a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tr->max_radius);
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type && !tr->M) {
107a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
108a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
109a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tr->M);CHKERRQ(ierr);
110a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(tr->M,tao->solution);CHKERRQ(ierr);
111a7e14dcfSSatish Balay   }
112a7e14dcfSSatish Balay 
113a7e14dcfSSatish Balay   /* Check convergence criteria */
114a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
115a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
11653506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
117a7e14dcfSSatish Balay   needH = 1;
118a7e14dcfSSatish Balay 
119a7e14dcfSSatish Balay   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
12053506e15SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay   /* Create vectors for the limited memory preconditioner */
123a7e14dcfSSatish Balay   if ((NTR_PC_BFGS == tr->pc_type) &&
124a7e14dcfSSatish Balay       (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) {
125a7e14dcfSSatish Balay     if (!tr->Diag) {
126a7e14dcfSSatish Balay         ierr = VecDuplicate(tao->solution, &tr->Diag);CHKERRQ(ierr);
127a7e14dcfSSatish Balay     }
128a7e14dcfSSatish Balay   }
129a7e14dcfSSatish Balay 
130a7e14dcfSSatish Balay   switch(tr->ksp_type) {
131a7e14dcfSSatish Balay   case NTR_KSP_NASH:
132a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
133a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
134a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
135a7e14dcfSSatish Balay     }
136a7e14dcfSSatish Balay     break;
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay   case NTR_KSP_STCG:
139a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
140a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
141a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
142a7e14dcfSSatish Balay     }
143a7e14dcfSSatish Balay     break;
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   default:
146a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
147a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
148a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
149a7e14dcfSSatish Balay     }
150a7e14dcfSSatish Balay     break;
151a7e14dcfSSatish Balay   }
152a7e14dcfSSatish Balay 
153a7e14dcfSSatish Balay   /*  Modify the preconditioner to use the bfgs approximation */
154a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
155a7e14dcfSSatish Balay   switch(tr->pc_type) {
156a7e14dcfSSatish Balay   case NTR_PC_NONE:
157a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
158a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
159a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
160a7e14dcfSSatish Balay     }
161a7e14dcfSSatish Balay     break;
162a7e14dcfSSatish Balay 
163a7e14dcfSSatish Balay   case NTR_PC_AHESS:
164a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
165a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
166a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
167a7e14dcfSSatish Balay     }
168a7e14dcfSSatish Balay     ierr = PCJacobiSetUseAbs(pc);CHKERRQ(ierr);
169a7e14dcfSSatish Balay     break;
170a7e14dcfSSatish Balay 
171a7e14dcfSSatish Balay   case NTR_PC_BFGS:
172a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
173a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
174a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
175a7e14dcfSSatish Balay     }
176a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
177a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, tr->M);CHKERRQ(ierr);
178a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
179a7e14dcfSSatish Balay     break;
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay   default:
182a7e14dcfSSatish Balay     /*  Use the pc method set by pc_type */
183a7e14dcfSSatish Balay     break;
184a7e14dcfSSatish Balay   }
185a7e14dcfSSatish Balay 
186a7e14dcfSSatish Balay   /*  Initialize trust-region radius */
187a7e14dcfSSatish Balay   switch(tr->init_type) {
188a7e14dcfSSatish Balay   case NTR_INIT_CONSTANT:
189a7e14dcfSSatish Balay     /*  Use the initial radius specified */
190a7e14dcfSSatish Balay     break;
191a7e14dcfSSatish Balay 
192a7e14dcfSSatish Balay   case NTR_INIT_INTERPOLATION:
193a7e14dcfSSatish Balay     /*  Use the initial radius specified */
194a7e14dcfSSatish Balay     max_radius = 0.0;
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
197a7e14dcfSSatish Balay       fmin = f;
198a7e14dcfSSatish Balay       sigma = 0.0;
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay       if (needH) {
201*ffad9901SBarry Smith         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
202a7e14dcfSSatish Balay         needH = 0;
203a7e14dcfSSatish Balay       }
204a7e14dcfSSatish Balay 
205a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
208a7e14dcfSSatish Balay         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
209a7e14dcfSSatish Balay         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
210a7e14dcfSSatish Balay 
211a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
212a7e14dcfSSatish Balay           tau = tr->gamma1_i;
213a7e14dcfSSatish Balay         }
214a7e14dcfSSatish Balay         else {
215a7e14dcfSSatish Balay           if (ftrial < fmin) {
216a7e14dcfSSatish Balay             fmin = ftrial;
217a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
218a7e14dcfSSatish Balay           }
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
221a7e14dcfSSatish Balay           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
224a7e14dcfSSatish Balay           actred = f - ftrial;
225a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
226a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tr->epsilon)) {
227a7e14dcfSSatish Balay             kappa = 1.0;
228a7e14dcfSSatish Balay           }
229a7e14dcfSSatish Balay           else {
230a7e14dcfSSatish Balay             kappa = actred / prered;
231a7e14dcfSSatish Balay           }
232a7e14dcfSSatish Balay 
233a7e14dcfSSatish Balay           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
234a7e14dcfSSatish Balay           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
235a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
236a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay           if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) {
239a7e14dcfSSatish Balay             /*  Great agreement */
240a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
241a7e14dcfSSatish Balay 
242a7e14dcfSSatish Balay             if (tau_max < 1.0) {
243a7e14dcfSSatish Balay               tau = tr->gamma3_i;
244a7e14dcfSSatish Balay             }
245a7e14dcfSSatish Balay             else if (tau_max > tr->gamma4_i) {
246a7e14dcfSSatish Balay               tau = tr->gamma4_i;
247a7e14dcfSSatish Balay             }
248a7e14dcfSSatish Balay             else {
249a7e14dcfSSatish Balay               tau = tau_max;
250a7e14dcfSSatish Balay             }
251a7e14dcfSSatish Balay           }
252a7e14dcfSSatish Balay           else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) {
253a7e14dcfSSatish Balay             /*  Good agreement */
254a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay             if (tau_max < tr->gamma2_i) {
257a7e14dcfSSatish Balay               tau = tr->gamma2_i;
258a7e14dcfSSatish Balay             }
259a7e14dcfSSatish Balay             else if (tau_max > tr->gamma3_i) {
260a7e14dcfSSatish Balay               tau = tr->gamma3_i;
261a7e14dcfSSatish Balay             }
262a7e14dcfSSatish Balay             else {
263a7e14dcfSSatish Balay               tau = tau_max;
264a7e14dcfSSatish Balay             }
265a7e14dcfSSatish Balay           }
266a7e14dcfSSatish Balay           else {
267a7e14dcfSSatish Balay             /*  Not good agreement */
268a7e14dcfSSatish Balay             if (tau_min > 1.0) {
269a7e14dcfSSatish Balay               tau = tr->gamma2_i;
270a7e14dcfSSatish Balay             }
271a7e14dcfSSatish Balay             else if (tau_max < tr->gamma1_i) {
272a7e14dcfSSatish Balay               tau = tr->gamma1_i;
273a7e14dcfSSatish Balay             }
274a7e14dcfSSatish Balay             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
275a7e14dcfSSatish Balay               tau = tr->gamma1_i;
276a7e14dcfSSatish Balay             }
277a7e14dcfSSatish Balay             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
278a7e14dcfSSatish Balay                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
279a7e14dcfSSatish Balay               tau = tau_1;
280a7e14dcfSSatish Balay             }
281a7e14dcfSSatish Balay             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
282a7e14dcfSSatish Balay                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
283a7e14dcfSSatish Balay               tau = tau_2;
284a7e14dcfSSatish Balay             }
285a7e14dcfSSatish Balay             else {
286a7e14dcfSSatish Balay               tau = tau_max;
287a7e14dcfSSatish Balay             }
288a7e14dcfSSatish Balay           }
289a7e14dcfSSatish Balay         }
290a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
291a7e14dcfSSatish Balay       }
292a7e14dcfSSatish Balay 
293a7e14dcfSSatish Balay       if (fmin < f) {
294a7e14dcfSSatish Balay         f = fmin;
295a7e14dcfSSatish Balay         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
296a7e14dcfSSatish Balay         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay         ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
299a7e14dcfSSatish Balay 
30053506e15SBarry Smith         if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
301a7e14dcfSSatish Balay         needH = 1;
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay         ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
304a7e14dcfSSatish Balay         if (reason != TAO_CONTINUE_ITERATING) {
305a7e14dcfSSatish Balay           PetscFunctionReturn(0);
306a7e14dcfSSatish Balay         }
307a7e14dcfSSatish Balay       }
308a7e14dcfSSatish Balay     }
309a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
310a7e14dcfSSatish Balay 
311a7e14dcfSSatish Balay     /*  Modify the radius if it is too large or small */
312a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tr->min_radius);
313a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
314a7e14dcfSSatish Balay     break;
315a7e14dcfSSatish Balay 
316a7e14dcfSSatish Balay   default:
317a7e14dcfSSatish Balay     /*  Norm of the first direction will initialize radius */
318a7e14dcfSSatish Balay     tao->trust = 0.0;
319a7e14dcfSSatish Balay     break;
320a7e14dcfSSatish Balay   }
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
323a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
324a7e14dcfSSatish Balay      since the function value may have decreased */
325a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type) {
326a7e14dcfSSatish Balay     if (f != 0.0) {
327a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
328a7e14dcfSSatish Balay     }
329a7e14dcfSSatish Balay     else {
330a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
331a7e14dcfSSatish Balay     }
332a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(tr->M,delta);CHKERRQ(ierr);
333a7e14dcfSSatish Balay   }
334a7e14dcfSSatish Balay 
335a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
336a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
337a7e14dcfSSatish Balay     ++iter;
338a7e14dcfSSatish Balay 
339a7e14dcfSSatish Balay     /* Compute the Hessian */
340a7e14dcfSSatish Balay     if (needH) {
341*ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
342a7e14dcfSSatish Balay       needH = 0;
343a7e14dcfSSatish Balay     }
344a7e14dcfSSatish Balay 
345a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type) {
346a7e14dcfSSatish Balay       if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) {
347a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner */
348a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, tr->Diag);CHKERRQ(ierr);
349a7e14dcfSSatish Balay         ierr = VecAbs(tr->Diag);CHKERRQ(ierr);
350a7e14dcfSSatish Balay         ierr = VecReciprocal(tr->Diag);CHKERRQ(ierr);
351a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(tr->M,tr->Diag);CHKERRQ(ierr);
352a7e14dcfSSatish Balay       }
353a7e14dcfSSatish Balay 
354a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
355a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
356a7e14dcfSSatish Balay       ++bfgsUpdates;
357a7e14dcfSSatish Balay     }
358a7e14dcfSSatish Balay 
359a7e14dcfSSatish Balay     while (reason == TAO_CONTINUE_ITERATING) {
36023ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
361a7e14dcfSSatish Balay 
362a7e14dcfSSatish Balay       /* Solve the trust region subproblem */
363a7e14dcfSSatish Balay       if (NTR_KSP_NASH == tr->ksp_type) {
364a7e14dcfSSatish Balay         ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
365a7e14dcfSSatish Balay         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
366a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
367a7e14dcfSSatish Balay         tao->ksp_its+=its;
368a7e14dcfSSatish Balay         ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
369a7e14dcfSSatish Balay       } else if (NTR_KSP_STCG == tr->ksp_type) {
370a7e14dcfSSatish Balay         ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
371a7e14dcfSSatish Balay         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
372a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
373a7e14dcfSSatish Balay         tao->ksp_its+=its;
374a7e14dcfSSatish Balay         ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
375a7e14dcfSSatish Balay       } else { /* NTR_KSP_GLTR */
376a7e14dcfSSatish Balay         ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
377a7e14dcfSSatish Balay         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
378a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
379a7e14dcfSSatish Balay         tao->ksp_its+=its;
380a7e14dcfSSatish Balay         ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
381a7e14dcfSSatish Balay       }
382a7e14dcfSSatish Balay 
383a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
384a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
385a7e14dcfSSatish Balay         if (norm_d > 0.0) {
386a7e14dcfSSatish Balay           tao->trust = norm_d;
387a7e14dcfSSatish Balay 
388a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
389a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
390a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
391a7e14dcfSSatish Balay         }
392a7e14dcfSSatish Balay         else {
393a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
394a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
395a7e14dcfSSatish Balay           tao->trust = tao->trust0;
396a7e14dcfSSatish Balay 
397a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
398a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
399a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
400a7e14dcfSSatish Balay 
401a7e14dcfSSatish Balay           if (NTR_KSP_NASH == tr->ksp_type) {
402a7e14dcfSSatish Balay             ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
403a7e14dcfSSatish Balay             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
404a7e14dcfSSatish Balay             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
405a7e14dcfSSatish Balay             tao->ksp_its+=its;
406a7e14dcfSSatish Balay             ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
407a7e14dcfSSatish Balay           } else if (NTR_KSP_STCG == tr->ksp_type) {
408a7e14dcfSSatish Balay             ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
409a7e14dcfSSatish Balay             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
410a7e14dcfSSatish Balay             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
411a7e14dcfSSatish Balay             tao->ksp_its+=its;
412a7e14dcfSSatish Balay             ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
413a7e14dcfSSatish Balay           } else { /* NTR_KSP_GLTR */
414a7e14dcfSSatish Balay             ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
415a7e14dcfSSatish Balay             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
416a7e14dcfSSatish Balay             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
417a7e14dcfSSatish Balay             tao->ksp_its+=its;
418a7e14dcfSSatish Balay             ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
419a7e14dcfSSatish Balay           }
420a7e14dcfSSatish Balay 
42153506e15SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
422a7e14dcfSSatish Balay         }
423a7e14dcfSSatish Balay       }
424a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
425a7e14dcfSSatish Balay       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
426a7e14dcfSSatish Balay       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&
427a7e14dcfSSatish Balay           (NTR_PC_BFGS == tr->pc_type) && (bfgsUpdates > 1)) {
428a7e14dcfSSatish Balay         /* Preconditioner is numerically indefinite; reset the
429a7e14dcfSSatish Balay            approximate if using BFGS preconditioning. */
430a7e14dcfSSatish Balay 
431a7e14dcfSSatish Balay         if (f != 0.0) {
432a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
433a7e14dcfSSatish Balay         }
434a7e14dcfSSatish Balay         else {
435a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
436a7e14dcfSSatish Balay         }
437a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(tr->M, delta);CHKERRQ(ierr);
438a7e14dcfSSatish Balay         ierr = MatLMVMReset(tr->M);CHKERRQ(ierr);
439a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
440a7e14dcfSSatish Balay         bfgsUpdates = 1;
441a7e14dcfSSatish Balay       }
442a7e14dcfSSatish Balay 
443a7e14dcfSSatish Balay       if (NTR_UPDATE_REDUCTION == tr->update_type) {
444a7e14dcfSSatish Balay         /* Get predicted reduction */
445a7e14dcfSSatish Balay         if (NTR_KSP_NASH == tr->ksp_type) {
446a7e14dcfSSatish Balay           ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
447a7e14dcfSSatish Balay         } else if (NTR_KSP_STCG == tr->ksp_type) {
448a7e14dcfSSatish Balay           ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
449a7e14dcfSSatish Balay         } else { /* gltr */
450a7e14dcfSSatish Balay           ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
451a7e14dcfSSatish Balay         }
452a7e14dcfSSatish Balay 
453a7e14dcfSSatish Balay         if (prered >= 0.0) {
454a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
455a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
456a7e14dcfSSatish Balay              be rejected! */
457a7e14dcfSSatish Balay           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
458a7e14dcfSSatish Balay         }
459a7e14dcfSSatish Balay         else {
460a7e14dcfSSatish Balay           /* Compute trial step and function value */
461a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
462a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
463a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
464a7e14dcfSSatish Balay 
465a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
466a7e14dcfSSatish Balay             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
467a7e14dcfSSatish Balay           } else {
468a7e14dcfSSatish Balay             /* Compute and actual reduction */
469a7e14dcfSSatish Balay             actred = f - ftrial;
470a7e14dcfSSatish Balay             prered = -prered;
471a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
472a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
473a7e14dcfSSatish Balay               kappa = 1.0;
474a7e14dcfSSatish Balay             }
475a7e14dcfSSatish Balay             else {
476a7e14dcfSSatish Balay               kappa = actred / prered;
477a7e14dcfSSatish Balay             }
478a7e14dcfSSatish Balay 
479a7e14dcfSSatish Balay             /* Accept or reject the step and update radius */
480a7e14dcfSSatish Balay             if (kappa < tr->eta1) {
481a7e14dcfSSatish Balay               /* Reject the step */
482a7e14dcfSSatish Balay               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
483a7e14dcfSSatish Balay             }
484a7e14dcfSSatish Balay             else {
485a7e14dcfSSatish Balay               /* Accept the step */
486a7e14dcfSSatish Balay               if (kappa < tr->eta2) {
487a7e14dcfSSatish Balay                 /* Marginal bad step */
488a7e14dcfSSatish Balay                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
489a7e14dcfSSatish Balay               }
490a7e14dcfSSatish Balay               else if (kappa < tr->eta3) {
491a7e14dcfSSatish Balay                 /* Reasonable step */
492a7e14dcfSSatish Balay                 tao->trust = tr->alpha3 * tao->trust;
493a7e14dcfSSatish Balay               }
494a7e14dcfSSatish Balay               else if (kappa < tr->eta4) {
495a7e14dcfSSatish Balay                 /* Good step */
496a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
497a7e14dcfSSatish Balay               }
498a7e14dcfSSatish Balay               else {
499a7e14dcfSSatish Balay                 /* Very good step */
500a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
501a7e14dcfSSatish Balay               }
502a7e14dcfSSatish Balay               break;
503a7e14dcfSSatish Balay             }
504a7e14dcfSSatish Balay           }
505a7e14dcfSSatish Balay         }
506a7e14dcfSSatish Balay       }
507a7e14dcfSSatish Balay       else {
508a7e14dcfSSatish Balay         /* Get predicted reduction */
509a7e14dcfSSatish Balay         if (NTR_KSP_NASH == tr->ksp_type) {
510a7e14dcfSSatish Balay           ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
511a7e14dcfSSatish Balay         } else if (NTR_KSP_STCG == tr->ksp_type) {
512a7e14dcfSSatish Balay           ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
513a7e14dcfSSatish Balay         } else { /* gltr */
514a7e14dcfSSatish Balay           ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
515a7e14dcfSSatish Balay         }
516a7e14dcfSSatish Balay 
517a7e14dcfSSatish Balay         if (prered >= 0.0) {
518a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
519a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
520a7e14dcfSSatish Balay              be rejected! */
521a7e14dcfSSatish Balay           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
522a7e14dcfSSatish Balay         }
523a7e14dcfSSatish Balay         else {
524a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
525a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
526a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
527a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
528a7e14dcfSSatish Balay             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
529a7e14dcfSSatish Balay           }
530a7e14dcfSSatish Balay           else {
531a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
532a7e14dcfSSatish Balay             actred = f - ftrial;
533a7e14dcfSSatish Balay             prered = -prered;
534a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
535a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
536a7e14dcfSSatish Balay               kappa = 1.0;
537a7e14dcfSSatish Balay             }
538a7e14dcfSSatish Balay             else {
539a7e14dcfSSatish Balay               kappa = actred / prered;
540a7e14dcfSSatish Balay             }
541a7e14dcfSSatish Balay 
542a7e14dcfSSatish Balay             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
543a7e14dcfSSatish Balay             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
544a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
545a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
546a7e14dcfSSatish Balay 
547a7e14dcfSSatish Balay             if (kappa >= 1.0 - tr->mu1) {
548a7e14dcfSSatish Balay               /* Great agreement; accept step and update radius */
549a7e14dcfSSatish Balay               if (tau_max < 1.0) {
550a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
551a7e14dcfSSatish Balay               }
552a7e14dcfSSatish Balay               else if (tau_max > tr->gamma4) {
553a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
554a7e14dcfSSatish Balay               }
555a7e14dcfSSatish Balay               else {
556a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
557a7e14dcfSSatish Balay               }
558a7e14dcfSSatish Balay               break;
559a7e14dcfSSatish Balay             }
560a7e14dcfSSatish Balay             else if (kappa >= 1.0 - tr->mu2) {
561a7e14dcfSSatish Balay               /* Good agreement */
562a7e14dcfSSatish Balay 
563a7e14dcfSSatish Balay               if (tau_max < tr->gamma2) {
564a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
565a7e14dcfSSatish Balay               }
566a7e14dcfSSatish Balay               else if (tau_max > tr->gamma3) {
567a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
568a7e14dcfSSatish Balay               }
569a7e14dcfSSatish Balay               else if (tau_max < 1.0) {
570a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
571a7e14dcfSSatish Balay               }
572a7e14dcfSSatish Balay               else {
573a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
574a7e14dcfSSatish Balay               }
575a7e14dcfSSatish Balay               break;
576a7e14dcfSSatish Balay             }
577a7e14dcfSSatish Balay             else {
578a7e14dcfSSatish Balay               /* Not good agreement */
579a7e14dcfSSatish Balay               if (tau_min > 1.0) {
580a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
581a7e14dcfSSatish Balay               }
582a7e14dcfSSatish Balay               else if (tau_max < tr->gamma1) {
583a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
584a7e14dcfSSatish Balay               }
585a7e14dcfSSatish Balay               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
586a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
587a7e14dcfSSatish Balay               }
588a7e14dcfSSatish Balay               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
589a7e14dcfSSatish Balay                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
590a7e14dcfSSatish Balay                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
591a7e14dcfSSatish Balay               }
592a7e14dcfSSatish Balay               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
593a7e14dcfSSatish Balay                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
594a7e14dcfSSatish Balay                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
595a7e14dcfSSatish Balay               }
596a7e14dcfSSatish Balay               else {
597a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
598a7e14dcfSSatish Balay               }
599a7e14dcfSSatish Balay             }
600a7e14dcfSSatish Balay           }
601a7e14dcfSSatish Balay         }
602a7e14dcfSSatish Balay       }
603a7e14dcfSSatish Balay 
604a7e14dcfSSatish Balay       /* The step computed was not good and the radius was decreased.
605a7e14dcfSSatish Balay          Monitor the radius to terminate. */
606a7e14dcfSSatish Balay       ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
607a7e14dcfSSatish Balay     }
608a7e14dcfSSatish Balay 
609a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
610a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
611a7e14dcfSSatish Balay 
612a7e14dcfSSatish Balay     if (reason == TAO_CONTINUE_ITERATING) {
613a7e14dcfSSatish Balay       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
614a7e14dcfSSatish Balay       f = ftrial;
615a7e14dcfSSatish Balay       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);
616a7e14dcfSSatish Balay       ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
61753506e15SBarry Smith       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
618a7e14dcfSSatish Balay       needH = 1;
619a7e14dcfSSatish Balay       ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
620a7e14dcfSSatish Balay     }
621a7e14dcfSSatish Balay   }
622a7e14dcfSSatish Balay   PetscFunctionReturn(0);
623a7e14dcfSSatish Balay }
624a7e14dcfSSatish Balay 
625a7e14dcfSSatish Balay /*------------------------------------------------------------*/
626a7e14dcfSSatish Balay #undef __FUNCT__
627a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTR"
628441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao)
629a7e14dcfSSatish Balay {
630a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
631a7e14dcfSSatish Balay   PetscErrorCode ierr;
632a7e14dcfSSatish Balay 
633a7e14dcfSSatish Balay   PetscFunctionBegin;
634a7e14dcfSSatish Balay 
635a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
636a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
637a7e14dcfSSatish Balay   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
638a7e14dcfSSatish Balay 
639a7e14dcfSSatish Balay   tr->Diag = 0;
640a7e14dcfSSatish Balay   tr->M = 0;
641a7e14dcfSSatish Balay 
642a7e14dcfSSatish Balay 
643a7e14dcfSSatish Balay   PetscFunctionReturn(0);
644a7e14dcfSSatish Balay }
645a7e14dcfSSatish Balay 
646a7e14dcfSSatish Balay /*------------------------------------------------------------*/
647a7e14dcfSSatish Balay #undef __FUNCT__
648a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTR"
649441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTR(Tao tao)
650a7e14dcfSSatish Balay {
651a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
652a7e14dcfSSatish Balay   PetscErrorCode ierr;
653a7e14dcfSSatish Balay 
654a7e14dcfSSatish Balay   PetscFunctionBegin;
655a7e14dcfSSatish Balay   if (tao->setupcalled) {
656a7e14dcfSSatish Balay     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
657a7e14dcfSSatish Balay   }
658a7e14dcfSSatish Balay   ierr = MatDestroy(&tr->M);CHKERRQ(ierr);
659a7e14dcfSSatish Balay   ierr = VecDestroy(&tr->Diag);CHKERRQ(ierr);
660a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
661a7e14dcfSSatish Balay   PetscFunctionReturn(0);
662a7e14dcfSSatish Balay }
663a7e14dcfSSatish Balay 
664a7e14dcfSSatish Balay /*------------------------------------------------------------*/
665a7e14dcfSSatish Balay #undef __FUNCT__
666a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTR"
667441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(Tao tao)
668a7e14dcfSSatish Balay {
669a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
670a7e14dcfSSatish Balay   PetscErrorCode ierr;
671a7e14dcfSSatish Balay 
672a7e14dcfSSatish Balay   PetscFunctionBegin;
673a7e14dcfSSatish Balay   ierr = PetscOptionsHead("Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
674a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_ntr_ksp_type", "ksp type", "", NTR_KSP, NTR_KSP_TYPES, NTR_KSP[tr->ksp_type], &tr->ksp_type, 0);CHKERRQ(ierr);
675a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_ntr_pc_type", "pc type", "", NTR_PC, NTR_PC_TYPES, NTR_PC[tr->pc_type], &tr->pc_type, 0);CHKERRQ(ierr);
676a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_ntr_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[tr->bfgs_scale_type], &tr->bfgs_scale_type, 0);CHKERRQ(ierr);
677a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_ntr_init_type", "tao->trust initialization type", "", NTR_INIT, NTR_INIT_TYPES, NTR_INIT[tr->init_type], &tr->init_type, 0);CHKERRQ(ierr);
678a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_ntr_update_type", "radius update type", "", NTR_UPDATE, NTR_UPDATE_TYPES, NTR_UPDATE[tr->update_type], &tr->update_type, 0);CHKERRQ(ierr);
679a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1, 0);CHKERRQ(ierr);
680a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2, 0);CHKERRQ(ierr);
681a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3, 0);CHKERRQ(ierr);
682a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4, 0);CHKERRQ(ierr);
683a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1, 0);CHKERRQ(ierr);
684a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2, 0);CHKERRQ(ierr);
685a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3, 0);CHKERRQ(ierr);
686a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4, 0);CHKERRQ(ierr);
687a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5, 0);CHKERRQ(ierr);
688a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1, 0);CHKERRQ(ierr);
689a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2, 0);CHKERRQ(ierr);
690a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1, 0);CHKERRQ(ierr);
691a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2, 0);CHKERRQ(ierr);
692a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3, 0);CHKERRQ(ierr);
693a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4, 0);CHKERRQ(ierr);
694a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta, 0);CHKERRQ(ierr);
695a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i, 0);CHKERRQ(ierr);
696a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i, 0);CHKERRQ(ierr);
697a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i, 0);CHKERRQ(ierr);
698a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i, 0);CHKERRQ(ierr);
699a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i, 0);CHKERRQ(ierr);
700a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i, 0);CHKERRQ(ierr);
701a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i, 0);CHKERRQ(ierr);
702a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius, 0);CHKERRQ(ierr);
703a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius, 0);CHKERRQ(ierr);
704a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon, 0);CHKERRQ(ierr);
705a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
706a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
707a7e14dcfSSatish Balay   PetscFunctionReturn(0);
708a7e14dcfSSatish Balay }
709a7e14dcfSSatish Balay 
710a7e14dcfSSatish Balay /*------------------------------------------------------------*/
711a7e14dcfSSatish Balay #undef __FUNCT__
712a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTR"
713441846f8SBarry Smith static PetscErrorCode TaoView_NTR(Tao tao, PetscViewer viewer)
714a7e14dcfSSatish Balay {
715a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
716a7e14dcfSSatish Balay   PetscErrorCode ierr;
717a7e14dcfSSatish Balay   PetscInt       nrejects;
718a7e14dcfSSatish Balay   PetscBool      isascii;
71953506e15SBarry Smith 
720a7e14dcfSSatish Balay   PetscFunctionBegin;
721a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
722a7e14dcfSSatish Balay   if (isascii) {
723a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
724a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type && tr->M) {
725a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(tr->M, &nrejects);CHKERRQ(ierr);
726a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr);
727a7e14dcfSSatish Balay     }
728a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
729a7e14dcfSSatish Balay   }
730a7e14dcfSSatish Balay   PetscFunctionReturn(0);
731a7e14dcfSSatish Balay }
732a7e14dcfSSatish Balay 
733a7e14dcfSSatish Balay /*------------------------------------------------------------*/
734a7e14dcfSSatish Balay EXTERN_C_BEGIN
735a7e14dcfSSatish Balay #undef __FUNCT__
736a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTR"
737441846f8SBarry Smith PetscErrorCode TaoCreate_NTR(Tao tao)
738a7e14dcfSSatish Balay {
739a7e14dcfSSatish Balay   TAO_NTR *tr;
740a7e14dcfSSatish Balay   PetscErrorCode ierr;
741a7e14dcfSSatish Balay 
742a7e14dcfSSatish Balay   PetscFunctionBegin;
743a7e14dcfSSatish Balay 
7443c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
745a7e14dcfSSatish Balay 
746a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTR;
747a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTR;
748a7e14dcfSSatish Balay   tao->ops->view = TaoView_NTR;
749a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
750a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTR;
751a7e14dcfSSatish Balay 
752a7e14dcfSSatish Balay   tao->max_it = 50;
7536f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7546f4723b1SBarry Smith   tao->fatol = 1e-5;
7556f4723b1SBarry Smith   tao->frtol = 1e-5;
7566f4723b1SBarry Smith #else
757a7e14dcfSSatish Balay   tao->fatol = 1e-10;
758a7e14dcfSSatish Balay   tao->frtol = 1e-10;
7596f4723b1SBarry Smith #endif
760a7e14dcfSSatish Balay   tao->data = (void*)tr;
761a7e14dcfSSatish Balay 
762a7e14dcfSSatish Balay   tao->trust0 = 100.0;
763a7e14dcfSSatish Balay 
764a7e14dcfSSatish Balay   /*  Standard trust region update parameters */
765a7e14dcfSSatish Balay   tr->eta1 = 1.0e-4;
766a7e14dcfSSatish Balay   tr->eta2 = 0.25;
767a7e14dcfSSatish Balay   tr->eta3 = 0.50;
768a7e14dcfSSatish Balay   tr->eta4 = 0.90;
769a7e14dcfSSatish Balay 
770a7e14dcfSSatish Balay   tr->alpha1 = 0.25;
771a7e14dcfSSatish Balay   tr->alpha2 = 0.50;
772a7e14dcfSSatish Balay   tr->alpha3 = 1.00;
773a7e14dcfSSatish Balay   tr->alpha4 = 2.00;
774a7e14dcfSSatish Balay   tr->alpha5 = 4.00;
775a7e14dcfSSatish Balay 
776a7e14dcfSSatish Balay   /*  Interpolation parameters */
777a7e14dcfSSatish Balay   tr->mu1_i = 0.35;
778a7e14dcfSSatish Balay   tr->mu2_i = 0.50;
779a7e14dcfSSatish Balay 
780a7e14dcfSSatish Balay   tr->gamma1_i = 0.0625;
781a7e14dcfSSatish Balay   tr->gamma2_i = 0.50;
782a7e14dcfSSatish Balay   tr->gamma3_i = 2.00;
783a7e14dcfSSatish Balay   tr->gamma4_i = 5.00;
784a7e14dcfSSatish Balay 
785a7e14dcfSSatish Balay   tr->theta_i = 0.25;
786a7e14dcfSSatish Balay 
787a7e14dcfSSatish Balay   /*  Interpolation trust region update parameters */
788a7e14dcfSSatish Balay   tr->mu1 = 0.10;
789a7e14dcfSSatish Balay   tr->mu2 = 0.50;
790a7e14dcfSSatish Balay 
791a7e14dcfSSatish Balay   tr->gamma1 = 0.25;
792a7e14dcfSSatish Balay   tr->gamma2 = 0.50;
793a7e14dcfSSatish Balay   tr->gamma3 = 2.00;
794a7e14dcfSSatish Balay   tr->gamma4 = 4.00;
795a7e14dcfSSatish Balay 
796a7e14dcfSSatish Balay   tr->theta = 0.05;
797a7e14dcfSSatish Balay 
798a7e14dcfSSatish Balay   tr->min_radius = 1.0e-10;
799a7e14dcfSSatish Balay   tr->max_radius = 1.0e10;
800a7e14dcfSSatish Balay   tr->epsilon = 1.0e-6;
801a7e14dcfSSatish Balay 
802a7e14dcfSSatish Balay   tr->ksp_type        = NTR_KSP_STCG;
803a7e14dcfSSatish Balay   tr->pc_type         = NTR_PC_BFGS;
804a7e14dcfSSatish Balay   tr->bfgs_scale_type = BFGS_SCALE_AHESS;
805a7e14dcfSSatish Balay   tr->init_type       = NTR_INIT_INTERPOLATION;
806a7e14dcfSSatish Balay   tr->update_type     = NTR_UPDATE_REDUCTION;
807a7e14dcfSSatish Balay 
808a7e14dcfSSatish Balay 
809a7e14dcfSSatish Balay   /* Set linear solver to default for trust region */
810a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
811a7e14dcfSSatish Balay 
812a7e14dcfSSatish Balay   PetscFunctionReturn(0);
813a7e14dcfSSatish Balay 
814a7e14dcfSSatish Balay 
815a7e14dcfSSatish Balay }
816a7e14dcfSSatish Balay EXTERN_C_END
817a7e14dcfSSatish Balay 
818a7e14dcfSSatish Balay 
819a7e14dcfSSatish Balay #undef __FUNCT__
820a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
821a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
822a7e14dcfSSatish Balay {
823a7e14dcfSSatish Balay     PetscErrorCode ierr;
824a7e14dcfSSatish Balay     Mat M;
825a7e14dcfSSatish Balay     PetscFunctionBegin;
826a7e14dcfSSatish Balay     PetscValidHeaderSpecific(pc,PC_CLASSID,1);
827a7e14dcfSSatish Balay     PetscValidHeaderSpecific(b,VEC_CLASSID,2);
828a7e14dcfSSatish Balay     PetscValidHeaderSpecific(x,VEC_CLASSID,3);
829a7e14dcfSSatish Balay     ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
830a7e14dcfSSatish Balay     ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
831a7e14dcfSSatish Balay     PetscFunctionReturn(0);
832a7e14dcfSSatish Balay }
833