xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision 5d5277661bc9d5de49c003e79a2b7124bf8a2eb4)
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>
6af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
7af0996ceSBarry 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           bfgsUpdates = 0;
87a7e14dcfSSatish Balay   PetscInt           needH;
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay   PetscInt           i_max = 5;
90a7e14dcfSSatish Balay   PetscInt           j_max = 1;
91a7e14dcfSSatish Balay   PetscInt           i, j, N, n, its;
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay   PetscFunctionBegin;
94a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
95a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr);
96a7e14dcfSSatish Balay   }
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   tao->trust = tao->trust0;
99a7e14dcfSSatish Balay 
100a7e14dcfSSatish Balay   /* Modify the radius if it is too large or small */
101a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tr->min_radius);
102a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tr->max_radius);
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type && !tr->M) {
106a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
107a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
108a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tr->M);CHKERRQ(ierr);
109a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(tr->M,tao->solution);CHKERRQ(ierr);
110a7e14dcfSSatish Balay   }
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay   /* Check convergence criteria */
113a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
114a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
11553506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
116a7e14dcfSSatish Balay   needH = 1;
117a7e14dcfSSatish Balay 
1188931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
11953506e15SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   /* Create vectors for the limited memory preconditioner */
122a7e14dcfSSatish Balay   if ((NTR_PC_BFGS == tr->pc_type) &&
123a7e14dcfSSatish Balay       (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) {
124a7e14dcfSSatish Balay     if (!tr->Diag) {
125a7e14dcfSSatish Balay         ierr = VecDuplicate(tao->solution, &tr->Diag);CHKERRQ(ierr);
126a7e14dcfSSatish Balay     }
127a7e14dcfSSatish Balay   }
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay   switch(tr->ksp_type) {
130a7e14dcfSSatish Balay   case NTR_KSP_NASH:
131a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
1321a1499c8SBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
133a7e14dcfSSatish Balay     break;
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay   case NTR_KSP_STCG:
136a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
1371a1499c8SBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
138a7e14dcfSSatish Balay     break;
139a7e14dcfSSatish Balay 
140a7e14dcfSSatish Balay   default:
141a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
1421a1499c8SBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
143a7e14dcfSSatish Balay     break;
144a7e14dcfSSatish Balay   }
145a7e14dcfSSatish Balay 
146a7e14dcfSSatish Balay   /*  Modify the preconditioner to use the bfgs approximation */
147a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
148a7e14dcfSSatish Balay   switch(tr->pc_type) {
149a7e14dcfSSatish Balay   case NTR_PC_NONE:
150a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1511a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
152a7e14dcfSSatish Balay     break;
153a7e14dcfSSatish Balay 
154a7e14dcfSSatish Balay   case NTR_PC_AHESS:
155a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1561a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
157baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
158a7e14dcfSSatish Balay     break;
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay   case NTR_PC_BFGS:
161a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
1621a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
163a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
164a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, tr->M);CHKERRQ(ierr);
165a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
166a7e14dcfSSatish Balay     break;
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   default:
169a7e14dcfSSatish Balay     /*  Use the pc method set by pc_type */
170a7e14dcfSSatish Balay     break;
171a7e14dcfSSatish Balay   }
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay   /*  Initialize trust-region radius */
174a7e14dcfSSatish Balay   switch(tr->init_type) {
175a7e14dcfSSatish Balay   case NTR_INIT_CONSTANT:
176a7e14dcfSSatish Balay     /*  Use the initial radius specified */
177a7e14dcfSSatish Balay     break;
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay   case NTR_INIT_INTERPOLATION:
180a7e14dcfSSatish Balay     /*  Use the initial radius specified */
181a7e14dcfSSatish Balay     max_radius = 0.0;
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
184a7e14dcfSSatish Balay       fmin = f;
185a7e14dcfSSatish Balay       sigma = 0.0;
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay       if (needH) {
188ffad9901SBarry Smith         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
189a7e14dcfSSatish Balay         needH = 0;
190a7e14dcfSSatish Balay       }
191a7e14dcfSSatish Balay 
192a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
195a7e14dcfSSatish Balay         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
196a7e14dcfSSatish Balay         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
197a7e14dcfSSatish Balay 
198a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
199a7e14dcfSSatish Balay           tau = tr->gamma1_i;
200a7e14dcfSSatish Balay         }
201a7e14dcfSSatish Balay         else {
202a7e14dcfSSatish Balay           if (ftrial < fmin) {
203a7e14dcfSSatish Balay             fmin = ftrial;
204a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
205a7e14dcfSSatish Balay           }
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
208a7e14dcfSSatish Balay           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
209a7e14dcfSSatish Balay 
210a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
211a7e14dcfSSatish Balay           actred = f - ftrial;
212a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
213a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tr->epsilon)) {
214a7e14dcfSSatish Balay             kappa = 1.0;
215a7e14dcfSSatish Balay           }
216a7e14dcfSSatish Balay           else {
217a7e14dcfSSatish Balay             kappa = actred / prered;
218a7e14dcfSSatish Balay           }
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
221a7e14dcfSSatish Balay           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
222a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
223a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay           if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) {
226a7e14dcfSSatish Balay             /*  Great agreement */
227a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
228a7e14dcfSSatish Balay 
229a7e14dcfSSatish Balay             if (tau_max < 1.0) {
230a7e14dcfSSatish Balay               tau = tr->gamma3_i;
231a7e14dcfSSatish Balay             }
232a7e14dcfSSatish Balay             else if (tau_max > tr->gamma4_i) {
233a7e14dcfSSatish Balay               tau = tr->gamma4_i;
234a7e14dcfSSatish Balay             }
235a7e14dcfSSatish Balay             else {
236a7e14dcfSSatish Balay               tau = tau_max;
237a7e14dcfSSatish Balay             }
238a7e14dcfSSatish Balay           }
239a7e14dcfSSatish Balay           else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) {
240a7e14dcfSSatish Balay             /*  Good agreement */
241a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay             if (tau_max < tr->gamma2_i) {
244a7e14dcfSSatish Balay               tau = tr->gamma2_i;
245a7e14dcfSSatish Balay             }
246a7e14dcfSSatish Balay             else if (tau_max > tr->gamma3_i) {
247a7e14dcfSSatish Balay               tau = tr->gamma3_i;
248a7e14dcfSSatish Balay             }
249a7e14dcfSSatish Balay             else {
250a7e14dcfSSatish Balay               tau = tau_max;
251a7e14dcfSSatish Balay             }
252a7e14dcfSSatish Balay           }
253a7e14dcfSSatish Balay           else {
254a7e14dcfSSatish Balay             /*  Not good agreement */
255a7e14dcfSSatish Balay             if (tau_min > 1.0) {
256a7e14dcfSSatish Balay               tau = tr->gamma2_i;
257a7e14dcfSSatish Balay             }
258a7e14dcfSSatish Balay             else if (tau_max < tr->gamma1_i) {
259a7e14dcfSSatish Balay               tau = tr->gamma1_i;
260a7e14dcfSSatish Balay             }
261a7e14dcfSSatish Balay             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
262a7e14dcfSSatish Balay               tau = tr->gamma1_i;
263a7e14dcfSSatish Balay             }
264a7e14dcfSSatish Balay             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
265a7e14dcfSSatish Balay                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
266a7e14dcfSSatish Balay               tau = tau_1;
267a7e14dcfSSatish Balay             }
268a7e14dcfSSatish Balay             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
269a7e14dcfSSatish Balay                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
270a7e14dcfSSatish Balay               tau = tau_2;
271a7e14dcfSSatish Balay             }
272a7e14dcfSSatish Balay             else {
273a7e14dcfSSatish Balay               tau = tau_max;
274a7e14dcfSSatish Balay             }
275a7e14dcfSSatish Balay           }
276a7e14dcfSSatish Balay         }
277a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
278a7e14dcfSSatish Balay       }
279a7e14dcfSSatish Balay 
280a7e14dcfSSatish Balay       if (fmin < f) {
281a7e14dcfSSatish Balay         f = fmin;
282a7e14dcfSSatish Balay         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
283a7e14dcfSSatish Balay         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
284a7e14dcfSSatish Balay 
285a7e14dcfSSatish Balay         ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
286a7e14dcfSSatish Balay 
28753506e15SBarry Smith         if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
288a7e14dcfSSatish Balay         needH = 1;
289a7e14dcfSSatish Balay 
2908931d482SJason Sarich         ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
291a7e14dcfSSatish Balay         if (reason != TAO_CONTINUE_ITERATING) {
292a7e14dcfSSatish Balay           PetscFunctionReturn(0);
293a7e14dcfSSatish Balay         }
294a7e14dcfSSatish Balay       }
295a7e14dcfSSatish Balay     }
296a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay     /*  Modify the radius if it is too large or small */
299a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tr->min_radius);
300a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
301a7e14dcfSSatish Balay     break;
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay   default:
304a7e14dcfSSatish Balay     /*  Norm of the first direction will initialize radius */
305a7e14dcfSSatish Balay     tao->trust = 0.0;
306a7e14dcfSSatish Balay     break;
307a7e14dcfSSatish Balay   }
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
310a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
311a7e14dcfSSatish Balay      since the function value may have decreased */
312a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type) {
313a7e14dcfSSatish Balay     if (f != 0.0) {
314a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
315a7e14dcfSSatish Balay     }
316a7e14dcfSSatish Balay     else {
317a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
318a7e14dcfSSatish Balay     }
319a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(tr->M,delta);CHKERRQ(ierr);
320a7e14dcfSSatish Balay   }
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
323a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
3248931d482SJason Sarich     ++tao->niter;
325ae93cb3cSJason Sarich     tao->ksp_its=0;
326a7e14dcfSSatish Balay     /* Compute the Hessian */
327a7e14dcfSSatish Balay     if (needH) {
328ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
329a7e14dcfSSatish Balay       needH = 0;
330a7e14dcfSSatish Balay     }
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type) {
333a7e14dcfSSatish Balay       if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) {
334a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner */
335a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, tr->Diag);CHKERRQ(ierr);
336a7e14dcfSSatish Balay         ierr = VecAbs(tr->Diag);CHKERRQ(ierr);
337a7e14dcfSSatish Balay         ierr = VecReciprocal(tr->Diag);CHKERRQ(ierr);
338a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(tr->M,tr->Diag);CHKERRQ(ierr);
339a7e14dcfSSatish Balay       }
340a7e14dcfSSatish Balay 
341a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
342a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
343a7e14dcfSSatish Balay       ++bfgsUpdates;
344a7e14dcfSSatish Balay     }
345a7e14dcfSSatish Balay 
346a7e14dcfSSatish Balay     while (reason == TAO_CONTINUE_ITERATING) {
34723ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
348a7e14dcfSSatish Balay 
349a7e14dcfSSatish Balay       /* Solve the trust region subproblem */
350a7e14dcfSSatish Balay       if (NTR_KSP_NASH == tr->ksp_type) {
351a7e14dcfSSatish Balay         ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
352a7e14dcfSSatish Balay         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
353a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
354a7e14dcfSSatish Balay         tao->ksp_its+=its;
355ae93cb3cSJason Sarich         tao->ksp_tot_its+=its;
356a7e14dcfSSatish Balay         ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
357a7e14dcfSSatish Balay       } else if (NTR_KSP_STCG == tr->ksp_type) {
358a7e14dcfSSatish Balay         ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
359a7e14dcfSSatish Balay         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
360a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
361a7e14dcfSSatish Balay         tao->ksp_its+=its;
362ae93cb3cSJason Sarich         tao->ksp_tot_its+=its;
363a7e14dcfSSatish Balay         ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
364a7e14dcfSSatish Balay       } else { /* NTR_KSP_GLTR */
365a7e14dcfSSatish Balay         ierr = KSPGLTRSetRadius(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;
3692d9aa51bSJason Sarich         tao->ksp_tot_its+=its;
370a7e14dcfSSatish Balay         ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
371a7e14dcfSSatish Balay       }
372a7e14dcfSSatish Balay 
373a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
374a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
375a7e14dcfSSatish Balay         if (norm_d > 0.0) {
376a7e14dcfSSatish Balay           tao->trust = norm_d;
377a7e14dcfSSatish Balay 
378a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
379a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
380a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
381a7e14dcfSSatish Balay         }
382a7e14dcfSSatish Balay         else {
383a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
384a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
385a7e14dcfSSatish Balay           tao->trust = tao->trust0;
386a7e14dcfSSatish Balay 
387a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
388a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
389a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
390a7e14dcfSSatish Balay 
391a7e14dcfSSatish Balay           if (NTR_KSP_NASH == tr->ksp_type) {
392a7e14dcfSSatish Balay             ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
393a7e14dcfSSatish Balay             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
394a7e14dcfSSatish Balay             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
395a7e14dcfSSatish Balay             tao->ksp_its+=its;
3962d9aa51bSJason Sarich             tao->ksp_tot_its+=its;
397a7e14dcfSSatish Balay             ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
398a7e14dcfSSatish Balay           } else if (NTR_KSP_STCG == tr->ksp_type) {
399a7e14dcfSSatish Balay             ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
400a7e14dcfSSatish Balay             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
401a7e14dcfSSatish Balay             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
402a7e14dcfSSatish Balay             tao->ksp_its+=its;
4032d9aa51bSJason Sarich             tao->ksp_tot_its+=its;
404a7e14dcfSSatish Balay             ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
405a7e14dcfSSatish Balay           } else { /* NTR_KSP_GLTR */
406a7e14dcfSSatish Balay             ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
407a7e14dcfSSatish Balay             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
408a7e14dcfSSatish Balay             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
409a7e14dcfSSatish Balay             tao->ksp_its+=its;
4102d9aa51bSJason Sarich             tao->ksp_tot_its+=its;
411a7e14dcfSSatish Balay             ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
412a7e14dcfSSatish Balay           }
413a7e14dcfSSatish Balay 
41453506e15SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
415a7e14dcfSSatish Balay         }
416a7e14dcfSSatish Balay       }
417a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
418a7e14dcfSSatish Balay       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
419a7e14dcfSSatish Balay       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&
420a7e14dcfSSatish Balay           (NTR_PC_BFGS == tr->pc_type) && (bfgsUpdates > 1)) {
421a7e14dcfSSatish Balay         /* Preconditioner is numerically indefinite; reset the
422a7e14dcfSSatish Balay            approximate if using BFGS preconditioning. */
423a7e14dcfSSatish Balay 
424a7e14dcfSSatish Balay         if (f != 0.0) {
425a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
426a7e14dcfSSatish Balay         }
427a7e14dcfSSatish Balay         else {
428a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
429a7e14dcfSSatish Balay         }
430a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(tr->M, delta);CHKERRQ(ierr);
431a7e14dcfSSatish Balay         ierr = MatLMVMReset(tr->M);CHKERRQ(ierr);
432a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
433a7e14dcfSSatish Balay         bfgsUpdates = 1;
434a7e14dcfSSatish Balay       }
435a7e14dcfSSatish Balay 
436a7e14dcfSSatish Balay       if (NTR_UPDATE_REDUCTION == tr->update_type) {
437a7e14dcfSSatish Balay         /* Get predicted reduction */
438a7e14dcfSSatish Balay         if (NTR_KSP_NASH == tr->ksp_type) {
439a7e14dcfSSatish Balay           ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
440a7e14dcfSSatish Balay         } else if (NTR_KSP_STCG == tr->ksp_type) {
441a7e14dcfSSatish Balay           ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
442a7e14dcfSSatish Balay         } else { /* gltr */
443a7e14dcfSSatish Balay           ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
444a7e14dcfSSatish Balay         }
445a7e14dcfSSatish Balay 
446a7e14dcfSSatish Balay         if (prered >= 0.0) {
447a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
448a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
449a7e14dcfSSatish Balay              be rejected! */
450a7e14dcfSSatish Balay           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
451a7e14dcfSSatish Balay         }
452a7e14dcfSSatish Balay         else {
453a7e14dcfSSatish Balay           /* Compute trial step and function value */
454a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
455a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
456a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
457a7e14dcfSSatish Balay 
458a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
459a7e14dcfSSatish Balay             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
460a7e14dcfSSatish Balay           } else {
461a7e14dcfSSatish Balay             /* Compute and actual reduction */
462a7e14dcfSSatish Balay             actred = f - ftrial;
463a7e14dcfSSatish Balay             prered = -prered;
464a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
465a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
466a7e14dcfSSatish Balay               kappa = 1.0;
467a7e14dcfSSatish Balay             }
468a7e14dcfSSatish Balay             else {
469a7e14dcfSSatish Balay               kappa = actred / prered;
470a7e14dcfSSatish Balay             }
471a7e14dcfSSatish Balay 
472a7e14dcfSSatish Balay             /* Accept or reject the step and update radius */
473a7e14dcfSSatish Balay             if (kappa < tr->eta1) {
474a7e14dcfSSatish Balay               /* Reject the step */
475a7e14dcfSSatish Balay               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
476a7e14dcfSSatish Balay             }
477a7e14dcfSSatish Balay             else {
478a7e14dcfSSatish Balay               /* Accept the step */
479a7e14dcfSSatish Balay               if (kappa < tr->eta2) {
480a7e14dcfSSatish Balay                 /* Marginal bad step */
481a7e14dcfSSatish Balay                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
482a7e14dcfSSatish Balay               }
483a7e14dcfSSatish Balay               else if (kappa < tr->eta3) {
484a7e14dcfSSatish Balay                 /* Reasonable step */
485a7e14dcfSSatish Balay                 tao->trust = tr->alpha3 * tao->trust;
486a7e14dcfSSatish Balay               }
487a7e14dcfSSatish Balay               else if (kappa < tr->eta4) {
488a7e14dcfSSatish Balay                 /* Good step */
489a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
490a7e14dcfSSatish Balay               }
491a7e14dcfSSatish Balay               else {
492a7e14dcfSSatish Balay                 /* Very good step */
493a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
494a7e14dcfSSatish Balay               }
495a7e14dcfSSatish Balay               break;
496a7e14dcfSSatish Balay             }
497a7e14dcfSSatish Balay           }
498a7e14dcfSSatish Balay         }
499a7e14dcfSSatish Balay       }
500a7e14dcfSSatish Balay       else {
501a7e14dcfSSatish Balay         /* Get predicted reduction */
502a7e14dcfSSatish Balay         if (NTR_KSP_NASH == tr->ksp_type) {
503a7e14dcfSSatish Balay           ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
504a7e14dcfSSatish Balay         } else if (NTR_KSP_STCG == tr->ksp_type) {
505a7e14dcfSSatish Balay           ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
506a7e14dcfSSatish Balay         } else { /* gltr */
507a7e14dcfSSatish Balay           ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
508a7e14dcfSSatish Balay         }
509a7e14dcfSSatish Balay 
510a7e14dcfSSatish Balay         if (prered >= 0.0) {
511a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
512a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
513a7e14dcfSSatish Balay              be rejected! */
514a7e14dcfSSatish Balay           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
515a7e14dcfSSatish Balay         }
516a7e14dcfSSatish Balay         else {
517a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
518a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
519a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
520a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
521a7e14dcfSSatish Balay             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
522a7e14dcfSSatish Balay           }
523a7e14dcfSSatish Balay           else {
524a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
525a7e14dcfSSatish Balay             actred = f - ftrial;
526a7e14dcfSSatish Balay             prered = -prered;
527a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
528a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
529a7e14dcfSSatish Balay               kappa = 1.0;
530a7e14dcfSSatish Balay             }
531a7e14dcfSSatish Balay             else {
532a7e14dcfSSatish Balay               kappa = actred / prered;
533a7e14dcfSSatish Balay             }
534a7e14dcfSSatish Balay 
535a7e14dcfSSatish Balay             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
536a7e14dcfSSatish Balay             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
537a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
538a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
539a7e14dcfSSatish Balay 
540a7e14dcfSSatish Balay             if (kappa >= 1.0 - tr->mu1) {
541a7e14dcfSSatish Balay               /* Great agreement; accept step and update radius */
542a7e14dcfSSatish Balay               if (tau_max < 1.0) {
543a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
544a7e14dcfSSatish Balay               }
545a7e14dcfSSatish Balay               else if (tau_max > tr->gamma4) {
546a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
547a7e14dcfSSatish Balay               }
548a7e14dcfSSatish Balay               else {
549a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
550a7e14dcfSSatish Balay               }
551a7e14dcfSSatish Balay               break;
552a7e14dcfSSatish Balay             }
553a7e14dcfSSatish Balay             else if (kappa >= 1.0 - tr->mu2) {
554a7e14dcfSSatish Balay               /* Good agreement */
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay               if (tau_max < tr->gamma2) {
557a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
558a7e14dcfSSatish Balay               }
559a7e14dcfSSatish Balay               else if (tau_max > tr->gamma3) {
560a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
561a7e14dcfSSatish Balay               }
562a7e14dcfSSatish Balay               else if (tau_max < 1.0) {
563a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
564a7e14dcfSSatish Balay               }
565a7e14dcfSSatish Balay               else {
566a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
567a7e14dcfSSatish Balay               }
568a7e14dcfSSatish Balay               break;
569a7e14dcfSSatish Balay             }
570a7e14dcfSSatish Balay             else {
571a7e14dcfSSatish Balay               /* Not good agreement */
572a7e14dcfSSatish Balay               if (tau_min > 1.0) {
573a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
574a7e14dcfSSatish Balay               }
575a7e14dcfSSatish Balay               else if (tau_max < tr->gamma1) {
576a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
577a7e14dcfSSatish Balay               }
578a7e14dcfSSatish Balay               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
579a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
580a7e14dcfSSatish Balay               }
581a7e14dcfSSatish Balay               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
582a7e14dcfSSatish Balay                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
583a7e14dcfSSatish Balay                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
584a7e14dcfSSatish Balay               }
585a7e14dcfSSatish Balay               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
586a7e14dcfSSatish Balay                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
587a7e14dcfSSatish Balay                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
588a7e14dcfSSatish Balay               }
589a7e14dcfSSatish Balay               else {
590a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
591a7e14dcfSSatish Balay               }
592a7e14dcfSSatish Balay             }
593a7e14dcfSSatish Balay           }
594a7e14dcfSSatish Balay         }
595a7e14dcfSSatish Balay       }
596a7e14dcfSSatish Balay 
597a7e14dcfSSatish Balay       /* The step computed was not good and the radius was decreased.
598a7e14dcfSSatish Balay          Monitor the radius to terminate. */
5998931d482SJason Sarich       ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
600a7e14dcfSSatish Balay     }
601a7e14dcfSSatish Balay 
602a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
603a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
604a7e14dcfSSatish Balay 
605a7e14dcfSSatish Balay     if (reason == TAO_CONTINUE_ITERATING) {
606a7e14dcfSSatish Balay       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
607a7e14dcfSSatish Balay       f = ftrial;
608302440fdSBarry Smith       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr);
609a7e14dcfSSatish Balay       ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
61053506e15SBarry Smith       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
611a7e14dcfSSatish Balay       needH = 1;
6128931d482SJason Sarich       ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
613a7e14dcfSSatish Balay     }
614a7e14dcfSSatish Balay   }
615a7e14dcfSSatish Balay   PetscFunctionReturn(0);
616a7e14dcfSSatish Balay }
617a7e14dcfSSatish Balay 
618a7e14dcfSSatish Balay /*------------------------------------------------------------*/
619a7e14dcfSSatish Balay #undef __FUNCT__
620a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTR"
621441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao)
622a7e14dcfSSatish Balay {
623a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
624a7e14dcfSSatish Balay   PetscErrorCode ierr;
625a7e14dcfSSatish Balay 
626a7e14dcfSSatish Balay   PetscFunctionBegin;
627a7e14dcfSSatish Balay 
628a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
629a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
630a7e14dcfSSatish Balay   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
631a7e14dcfSSatish Balay 
632a7e14dcfSSatish Balay   tr->Diag = 0;
633a7e14dcfSSatish Balay   tr->M = 0;
634a7e14dcfSSatish Balay 
635a7e14dcfSSatish Balay 
636a7e14dcfSSatish Balay   PetscFunctionReturn(0);
637a7e14dcfSSatish Balay }
638a7e14dcfSSatish Balay 
639a7e14dcfSSatish Balay /*------------------------------------------------------------*/
640a7e14dcfSSatish Balay #undef __FUNCT__
641a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTR"
642441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTR(Tao tao)
643a7e14dcfSSatish Balay {
644a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
645a7e14dcfSSatish Balay   PetscErrorCode ierr;
646a7e14dcfSSatish Balay 
647a7e14dcfSSatish Balay   PetscFunctionBegin;
648a7e14dcfSSatish Balay   if (tao->setupcalled) {
649a7e14dcfSSatish Balay     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
650a7e14dcfSSatish Balay   }
651a7e14dcfSSatish Balay   ierr = MatDestroy(&tr->M);CHKERRQ(ierr);
652a7e14dcfSSatish Balay   ierr = VecDestroy(&tr->Diag);CHKERRQ(ierr);
653a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
654a7e14dcfSSatish Balay   PetscFunctionReturn(0);
655a7e14dcfSSatish Balay }
656a7e14dcfSSatish Balay 
657a7e14dcfSSatish Balay /*------------------------------------------------------------*/
658a7e14dcfSSatish Balay #undef __FUNCT__
659a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTR"
6608c34d3f5SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(PetscOptions *PetscOptionsObject,Tao tao)
661a7e14dcfSSatish Balay {
662a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
663a7e14dcfSSatish Balay   PetscErrorCode ierr;
664a7e14dcfSSatish Balay 
665a7e14dcfSSatish Balay   PetscFunctionBegin;
6661a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
66794ae4db5SBarry Smith   ierr = PetscOptionsEList("-tao_ntr_ksp_type", "ksp type", "", NTR_KSP, NTR_KSP_TYPES, NTR_KSP[tr->ksp_type], &tr->ksp_type,NULL);CHKERRQ(ierr);
66894ae4db5SBarry Smith   ierr = PetscOptionsEList("-tao_ntr_pc_type", "pc type", "", NTR_PC, NTR_PC_TYPES, NTR_PC[tr->pc_type], &tr->pc_type,NULL);CHKERRQ(ierr);
66994ae4db5SBarry Smith   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,NULL);CHKERRQ(ierr);
67094ae4db5SBarry Smith   ierr = PetscOptionsEList("-tao_ntr_init_type", "tao->trust initialization type", "", NTR_INIT, NTR_INIT_TYPES, NTR_INIT[tr->init_type], &tr->init_type,NULL);CHKERRQ(ierr);
67194ae4db5SBarry Smith   ierr = PetscOptionsEList("-tao_ntr_update_type", "radius update type", "", NTR_UPDATE, NTR_UPDATE_TYPES, NTR_UPDATE[tr->update_type], &tr->update_type,NULL);CHKERRQ(ierr);
67294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr);
67394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr);
67494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr);
67594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr);
67694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr);
67794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr);
67894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr);
67994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr);
68094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr);
68194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr);
68294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr);
68394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr);
68494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr);
68594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr);
68694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr);
68794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr);
68894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr);
68994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr);
69094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr);
69194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr);
69294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr);
69394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr);
69494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr);
69594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr);
69694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr);
69794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr);
698a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
699a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
700a7e14dcfSSatish Balay   PetscFunctionReturn(0);
701a7e14dcfSSatish Balay }
702a7e14dcfSSatish Balay 
703a7e14dcfSSatish Balay /*------------------------------------------------------------*/
704a7e14dcfSSatish Balay #undef __FUNCT__
705a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTR"
706441846f8SBarry Smith static PetscErrorCode TaoView_NTR(Tao tao, PetscViewer viewer)
707a7e14dcfSSatish Balay {
708a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
709a7e14dcfSSatish Balay   PetscErrorCode ierr;
710a7e14dcfSSatish Balay   PetscInt       nrejects;
711a7e14dcfSSatish Balay   PetscBool      isascii;
71253506e15SBarry Smith 
713a7e14dcfSSatish Balay   PetscFunctionBegin;
714a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
715a7e14dcfSSatish Balay   if (isascii) {
716a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
717a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type && tr->M) {
718a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(tr->M, &nrejects);CHKERRQ(ierr);
719a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr);
720a7e14dcfSSatish Balay     }
721a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
722a7e14dcfSSatish Balay   }
723a7e14dcfSSatish Balay   PetscFunctionReturn(0);
724a7e14dcfSSatish Balay }
725a7e14dcfSSatish Balay 
726a7e14dcfSSatish Balay /*------------------------------------------------------------*/
7271522df2eSJason Sarich /*MC
7281522df2eSJason Sarich   TAONTR - Newton's method with trust region for unconstrained minimization.
7291522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system.
7301522df2eSJason Sarich   NTR expects a KSP solver with a trust region radius.
7311522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
7321522df2eSJason Sarich 
7331522df2eSJason Sarich   Options Database Keys:
7341522df2eSJason Sarich + -tao_ntr_ksp_type - "nash","stcg","gltr"
7351522df2eSJason Sarich . -tao_ntr_pc_type - "none","ahess","bfgs","petsc"
7361522df2eSJason Sarich . -tao_ntr_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs"
7371522df2eSJason Sarich . -tao_ntr_init_type - "constant","direction","interpolation"
7381522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation"
7391522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius
7401522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius
7411522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
7421522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor
7431522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor
7441522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor
7451522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor
7461522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor
7471522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor
7481522df2eSJason Sarich . -tao_ntr_theta_i - thetha1 interpolation init factor
7491522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor
7501522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor
7511522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor
7521522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor
7531522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor
7541522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor
7551522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor
7561522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
7571522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
7581522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update
7591522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update
7601522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update
7611522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update
7621522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update
7631522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update
7641522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update
7651522df2eSJason Sarich 
7661eb8069cSJason Sarich   Level: beginner
7671522df2eSJason Sarich M*/
7681522df2eSJason Sarich 
769a7e14dcfSSatish Balay #undef __FUNCT__
770a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTR"
771728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
772a7e14dcfSSatish Balay {
773a7e14dcfSSatish Balay   TAO_NTR *tr;
774a7e14dcfSSatish Balay   PetscErrorCode ierr;
775a7e14dcfSSatish Balay 
776a7e14dcfSSatish Balay   PetscFunctionBegin;
777a7e14dcfSSatish Balay 
7783c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
779a7e14dcfSSatish Balay 
780a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTR;
781a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTR;
782a7e14dcfSSatish Balay   tao->ops->view = TaoView_NTR;
783a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
784a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTR;
785a7e14dcfSSatish Balay 
786a7e14dcfSSatish Balay   tao->max_it = 50;
7876f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7886f4723b1SBarry Smith   tao->fatol = 1e-5;
7896f4723b1SBarry Smith   tao->frtol = 1e-5;
7906f4723b1SBarry Smith #else
791a7e14dcfSSatish Balay   tao->fatol = 1e-10;
792a7e14dcfSSatish Balay   tao->frtol = 1e-10;
7936f4723b1SBarry Smith #endif
794a7e14dcfSSatish Balay   tao->data = (void*)tr;
795a7e14dcfSSatish Balay 
796a7e14dcfSSatish Balay   tao->trust0 = 100.0;
797a7e14dcfSSatish Balay 
798a7e14dcfSSatish Balay   /*  Standard trust region update parameters */
799a7e14dcfSSatish Balay   tr->eta1 = 1.0e-4;
800a7e14dcfSSatish Balay   tr->eta2 = 0.25;
801a7e14dcfSSatish Balay   tr->eta3 = 0.50;
802a7e14dcfSSatish Balay   tr->eta4 = 0.90;
803a7e14dcfSSatish Balay 
804a7e14dcfSSatish Balay   tr->alpha1 = 0.25;
805a7e14dcfSSatish Balay   tr->alpha2 = 0.50;
806a7e14dcfSSatish Balay   tr->alpha3 = 1.00;
807a7e14dcfSSatish Balay   tr->alpha4 = 2.00;
808a7e14dcfSSatish Balay   tr->alpha5 = 4.00;
809a7e14dcfSSatish Balay 
810a7e14dcfSSatish Balay   /*  Interpolation parameters */
811a7e14dcfSSatish Balay   tr->mu1_i = 0.35;
812a7e14dcfSSatish Balay   tr->mu2_i = 0.50;
813a7e14dcfSSatish Balay 
814a7e14dcfSSatish Balay   tr->gamma1_i = 0.0625;
815a7e14dcfSSatish Balay   tr->gamma2_i = 0.50;
816a7e14dcfSSatish Balay   tr->gamma3_i = 2.00;
817a7e14dcfSSatish Balay   tr->gamma4_i = 5.00;
818a7e14dcfSSatish Balay 
819a7e14dcfSSatish Balay   tr->theta_i = 0.25;
820a7e14dcfSSatish Balay 
821a7e14dcfSSatish Balay   /*  Interpolation trust region update parameters */
822a7e14dcfSSatish Balay   tr->mu1 = 0.10;
823a7e14dcfSSatish Balay   tr->mu2 = 0.50;
824a7e14dcfSSatish Balay 
825a7e14dcfSSatish Balay   tr->gamma1 = 0.25;
826a7e14dcfSSatish Balay   tr->gamma2 = 0.50;
827a7e14dcfSSatish Balay   tr->gamma3 = 2.00;
828a7e14dcfSSatish Balay   tr->gamma4 = 4.00;
829a7e14dcfSSatish Balay 
830a7e14dcfSSatish Balay   tr->theta = 0.05;
831a7e14dcfSSatish Balay 
832a7e14dcfSSatish Balay   tr->min_radius = 1.0e-10;
833a7e14dcfSSatish Balay   tr->max_radius = 1.0e10;
834a7e14dcfSSatish Balay   tr->epsilon = 1.0e-6;
835a7e14dcfSSatish Balay 
836a7e14dcfSSatish Balay   tr->ksp_type        = NTR_KSP_STCG;
837a7e14dcfSSatish Balay   tr->pc_type         = NTR_PC_BFGS;
838a7e14dcfSSatish Balay   tr->bfgs_scale_type = BFGS_SCALE_AHESS;
839a7e14dcfSSatish Balay   tr->init_type       = NTR_INIT_INTERPOLATION;
840a7e14dcfSSatish Balay   tr->update_type     = NTR_UPDATE_REDUCTION;
841a7e14dcfSSatish Balay 
842a7e14dcfSSatish Balay 
843a7e14dcfSSatish Balay   /* Set linear solver to default for trust region */
844a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
845*5d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
846a7e14dcfSSatish Balay   PetscFunctionReturn(0);
847a7e14dcfSSatish Balay }
848a7e14dcfSSatish Balay 
849a7e14dcfSSatish Balay 
850a7e14dcfSSatish Balay #undef __FUNCT__
851a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
852a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
853a7e14dcfSSatish Balay {
854a7e14dcfSSatish Balay     PetscErrorCode ierr;
855a7e14dcfSSatish Balay     Mat M;
856a7e14dcfSSatish Balay     PetscFunctionBegin;
857a7e14dcfSSatish Balay     PetscValidHeaderSpecific(pc,PC_CLASSID,1);
858a7e14dcfSSatish Balay     PetscValidHeaderSpecific(b,VEC_CLASSID,2);
859a7e14dcfSSatish Balay     PetscValidHeaderSpecific(x,VEC_CLASSID,3);
860a7e14dcfSSatish Balay     ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
861a7e14dcfSSatish Balay     ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
862a7e14dcfSSatish Balay     PetscFunctionReturn(0);
863a7e14dcfSSatish Balay }
864