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