xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision ba7fe8fb864bc2e28543a5593f37aff2eb73db37)
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>
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay #define NTR_KSP_NASH    0
7a7e14dcfSSatish Balay #define NTR_KSP_STCG    1
8a7e14dcfSSatish Balay #define NTR_KSP_GLTR    2
9a7e14dcfSSatish Balay #define NTR_KSP_TYPES   3
10a7e14dcfSSatish Balay 
11a7e14dcfSSatish Balay #define NTR_PC_NONE     0
12a7e14dcfSSatish Balay #define NTR_PC_AHESS    1
13a7e14dcfSSatish Balay #define NTR_PC_BFGS     2
14a7e14dcfSSatish Balay #define NTR_PC_PETSC    3
15a7e14dcfSSatish Balay #define NTR_PC_TYPES    4
16a7e14dcfSSatish Balay 
17a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS   0
18a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS    1
19a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES   2
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay #define NTR_INIT_CONSTANT         0
22a7e14dcfSSatish Balay #define NTR_INIT_DIRECTION        1
23a7e14dcfSSatish Balay #define NTR_INIT_INTERPOLATION    2
24a7e14dcfSSatish Balay #define NTR_INIT_TYPES            3
25a7e14dcfSSatish Balay 
26a7e14dcfSSatish Balay #define NTR_UPDATE_REDUCTION      0
27a7e14dcfSSatish Balay #define NTR_UPDATE_INTERPOLATION  1
28a7e14dcfSSatish Balay #define NTR_UPDATE_TYPES          2
29a7e14dcfSSatish Balay 
3053506e15SBarry Smith static const char *NTR_KSP[64] = {  "nash", "stcg", "gltr"};
31a7e14dcfSSatish Balay 
3253506e15SBarry Smith static const char *NTR_PC[64] = {  "none", "ahess", "bfgs", "petsc"};
33a7e14dcfSSatish Balay 
3453506e15SBarry Smith static const char *BFGS_SCALE[64] = {  "ahess", "bfgs"};
35a7e14dcfSSatish Balay 
3653506e15SBarry Smith static const char *NTR_INIT[64] = {  "constant", "direction", "interpolation"};
37a7e14dcfSSatish Balay 
3853506e15SBarry Smith static const char *NTR_UPDATE[64] = {  "reduction", "interpolation"};
39a7e14dcfSSatish Balay 
40a7e14dcfSSatish Balay /*  Routine for BFGS preconditioner */
41a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec xin, Vec xout);
42a7e14dcfSSatish Balay 
43a7e14dcfSSatish Balay /*
44a7e14dcfSSatish Balay    TaoSolve_NTR - Implements Newton's Method with a trust region approach
45a7e14dcfSSatish Balay    for solving unconstrained minimization problems.
46a7e14dcfSSatish Balay 
47a7e14dcfSSatish Balay    The basic algorithm is taken from MINPACK-2 (dstrn).
48a7e14dcfSSatish Balay 
49a7e14dcfSSatish Balay    TaoSolve_NTR computes a local minimizer of a twice differentiable function
50a7e14dcfSSatish Balay    f by applying a trust region variant of Newton's method.  At each stage
51a7e14dcfSSatish Balay    of the algorithm, we use the prconditioned conjugate gradient method to
52a7e14dcfSSatish Balay    determine an approximate minimizer of the quadratic equation
53a7e14dcfSSatish Balay 
54a7e14dcfSSatish Balay         q(s) = <s, Hs + g>
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay    subject to the trust region constraint
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay         || s ||_M <= radius,
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay    where radius is the trust region radius and M is a symmetric positive
61a7e14dcfSSatish Balay    definite matrix (the preconditioner).  Here g is the gradient and H
62a7e14dcfSSatish Balay    is the Hessian matrix.
63a7e14dcfSSatish Balay 
64*ba7fe8fbSTodd Munson    Note:  TaoSolve_NTR MUST use the iterative solver KSPCGNASH, KSPCGSTCG,
65*ba7fe8fbSTodd Munson           or KSPCGGLTR.  Thus, we set KSPCGNASH, KSPCGSTCG, or KSPCGGLTR in this
66a7e14dcfSSatish Balay           routine regardless of what the user may have previously specified.
67a7e14dcfSSatish Balay */
68a7e14dcfSSatish Balay #undef __FUNCT__
69a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NTR"
70441846f8SBarry Smith static PetscErrorCode TaoSolve_NTR(Tao tao)
71a7e14dcfSSatish Balay {
72a7e14dcfSSatish Balay   TAO_NTR            *tr = (TAO_NTR *)tao->data;
73a7e14dcfSSatish Balay   PC                 pc;
74a7e14dcfSSatish Balay   KSPConvergedReason ksp_reason;
75e4cb33bbSBarry Smith   TaoConvergedReason reason;
76a7e14dcfSSatish Balay   PetscReal          fmin, ftrial, prered, actred, kappa, sigma, beta;
77a7e14dcfSSatish Balay   PetscReal          tau, tau_1, tau_2, tau_max, tau_min, max_radius;
78a7e14dcfSSatish Balay   PetscReal          f, gnorm;
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay   PetscReal          delta;
81a7e14dcfSSatish Balay   PetscReal          norm_d;
82a7e14dcfSSatish Balay   PetscErrorCode     ierr;
83a7e14dcfSSatish Balay   PetscInt           bfgsUpdates = 0;
84a7e14dcfSSatish Balay   PetscInt           needH;
85a7e14dcfSSatish Balay 
86a7e14dcfSSatish Balay   PetscInt           i_max = 5;
87a7e14dcfSSatish Balay   PetscInt           j_max = 1;
88a7e14dcfSSatish Balay   PetscInt           i, j, N, n, its;
89a7e14dcfSSatish Balay 
90a7e14dcfSSatish Balay   PetscFunctionBegin;
91a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
92a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr);
93a7e14dcfSSatish Balay   }
94a7e14dcfSSatish Balay 
95a7e14dcfSSatish Balay   tao->trust = tao->trust0;
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay   /* Modify the radius if it is too large or small */
98a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tr->min_radius);
99a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tr->max_radius);
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type && !tr->M) {
103a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
104a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
105a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tr->M);CHKERRQ(ierr);
106a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(tr->M,tao->solution);CHKERRQ(ierr);
107a7e14dcfSSatish Balay   }
108a7e14dcfSSatish Balay 
109a7e14dcfSSatish Balay   /* Check convergence criteria */
110a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
111a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
11253506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
113a7e14dcfSSatish Balay   needH = 1;
114a7e14dcfSSatish Balay 
1158931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
11653506e15SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   /* Create vectors for the limited memory preconditioner */
119a7e14dcfSSatish Balay   if ((NTR_PC_BFGS == tr->pc_type) &&
120a7e14dcfSSatish Balay       (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) {
121a7e14dcfSSatish Balay     if (!tr->Diag) {
122a7e14dcfSSatish Balay         ierr = VecDuplicate(tao->solution, &tr->Diag);CHKERRQ(ierr);
123a7e14dcfSSatish Balay     }
124a7e14dcfSSatish Balay   }
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay   switch(tr->ksp_type) {
127a7e14dcfSSatish Balay   case NTR_KSP_NASH:
128*ba7fe8fbSTodd Munson     ierr = KSPSetType(tao->ksp, KSPCGNASH);CHKERRQ(ierr);
1291a1499c8SBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
130a7e14dcfSSatish Balay     break;
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay   case NTR_KSP_STCG:
133*ba7fe8fbSTodd Munson     ierr = KSPSetType(tao->ksp, KSPCGSTCG);CHKERRQ(ierr);
1341a1499c8SBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
135a7e14dcfSSatish Balay     break;
136a7e14dcfSSatish Balay 
137a7e14dcfSSatish Balay   default:
138*ba7fe8fbSTodd Munson     ierr = KSPSetType(tao->ksp, KSPCGGLTR);CHKERRQ(ierr);
1391a1499c8SBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
140a7e14dcfSSatish Balay     break;
141a7e14dcfSSatish Balay   }
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay   /*  Modify the preconditioner to use the bfgs approximation */
144a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
145a7e14dcfSSatish Balay   switch(tr->pc_type) {
146a7e14dcfSSatish Balay   case NTR_PC_NONE:
147a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1481a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
149a7e14dcfSSatish Balay     break;
150a7e14dcfSSatish Balay 
151a7e14dcfSSatish Balay   case NTR_PC_AHESS:
152a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1531a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
154baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
155a7e14dcfSSatish Balay     break;
156a7e14dcfSSatish Balay 
157a7e14dcfSSatish Balay   case NTR_PC_BFGS:
158a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
1591a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
160a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
161a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, tr->M);CHKERRQ(ierr);
162a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
163a7e14dcfSSatish Balay     break;
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay   default:
166a7e14dcfSSatish Balay     /*  Use the pc method set by pc_type */
167a7e14dcfSSatish Balay     break;
168a7e14dcfSSatish Balay   }
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay   /*  Initialize trust-region radius */
171a7e14dcfSSatish Balay   switch(tr->init_type) {
172a7e14dcfSSatish Balay   case NTR_INIT_CONSTANT:
173a7e14dcfSSatish Balay     /*  Use the initial radius specified */
174a7e14dcfSSatish Balay     break;
175a7e14dcfSSatish Balay 
176a7e14dcfSSatish Balay   case NTR_INIT_INTERPOLATION:
177a7e14dcfSSatish Balay     /*  Use the initial radius specified */
178a7e14dcfSSatish Balay     max_radius = 0.0;
179a7e14dcfSSatish Balay 
180a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
181a7e14dcfSSatish Balay       fmin = f;
182a7e14dcfSSatish Balay       sigma = 0.0;
183a7e14dcfSSatish Balay 
184a7e14dcfSSatish Balay       if (needH) {
185ffad9901SBarry Smith         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
186a7e14dcfSSatish Balay         needH = 0;
187a7e14dcfSSatish Balay       }
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
192a7e14dcfSSatish Balay         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
193a7e14dcfSSatish Balay         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
196a7e14dcfSSatish Balay           tau = tr->gamma1_i;
197a7e14dcfSSatish Balay         }
198a7e14dcfSSatish Balay         else {
199a7e14dcfSSatish Balay           if (ftrial < fmin) {
200a7e14dcfSSatish Balay             fmin = ftrial;
201a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
202a7e14dcfSSatish Balay           }
203a7e14dcfSSatish Balay 
204a7e14dcfSSatish Balay           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
205a7e14dcfSSatish Balay           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
208a7e14dcfSSatish Balay           actred = f - ftrial;
209a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
210a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tr->epsilon)) {
211a7e14dcfSSatish Balay             kappa = 1.0;
212a7e14dcfSSatish Balay           }
213a7e14dcfSSatish Balay           else {
214a7e14dcfSSatish Balay             kappa = actred / prered;
215a7e14dcfSSatish Balay           }
216a7e14dcfSSatish Balay 
217a7e14dcfSSatish Balay           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
218a7e14dcfSSatish Balay           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
219a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
220a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
221a7e14dcfSSatish Balay 
222a7e14dcfSSatish Balay           if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) {
223a7e14dcfSSatish Balay             /*  Great agreement */
224a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay             if (tau_max < 1.0) {
227a7e14dcfSSatish Balay               tau = tr->gamma3_i;
228a7e14dcfSSatish Balay             }
229a7e14dcfSSatish Balay             else if (tau_max > tr->gamma4_i) {
230a7e14dcfSSatish Balay               tau = tr->gamma4_i;
231a7e14dcfSSatish Balay             }
232a7e14dcfSSatish Balay             else {
233a7e14dcfSSatish Balay               tau = tau_max;
234a7e14dcfSSatish Balay             }
235a7e14dcfSSatish Balay           }
236a7e14dcfSSatish Balay           else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) {
237a7e14dcfSSatish Balay             /*  Good agreement */
238a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay             if (tau_max < tr->gamma2_i) {
241a7e14dcfSSatish Balay               tau = tr->gamma2_i;
242a7e14dcfSSatish Balay             }
243a7e14dcfSSatish Balay             else if (tau_max > tr->gamma3_i) {
244a7e14dcfSSatish Balay               tau = tr->gamma3_i;
245a7e14dcfSSatish Balay             }
246a7e14dcfSSatish Balay             else {
247a7e14dcfSSatish Balay               tau = tau_max;
248a7e14dcfSSatish Balay             }
249a7e14dcfSSatish Balay           }
250a7e14dcfSSatish Balay           else {
251a7e14dcfSSatish Balay             /*  Not good agreement */
252a7e14dcfSSatish Balay             if (tau_min > 1.0) {
253a7e14dcfSSatish Balay               tau = tr->gamma2_i;
254a7e14dcfSSatish Balay             }
255a7e14dcfSSatish Balay             else if (tau_max < tr->gamma1_i) {
256a7e14dcfSSatish Balay               tau = tr->gamma1_i;
257a7e14dcfSSatish Balay             }
258a7e14dcfSSatish Balay             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
259a7e14dcfSSatish Balay               tau = tr->gamma1_i;
260a7e14dcfSSatish Balay             }
261a7e14dcfSSatish Balay             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
262a7e14dcfSSatish Balay                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
263a7e14dcfSSatish Balay               tau = tau_1;
264a7e14dcfSSatish Balay             }
265a7e14dcfSSatish Balay             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
266a7e14dcfSSatish Balay                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
267a7e14dcfSSatish Balay               tau = tau_2;
268a7e14dcfSSatish Balay             }
269a7e14dcfSSatish Balay             else {
270a7e14dcfSSatish Balay               tau = tau_max;
271a7e14dcfSSatish Balay             }
272a7e14dcfSSatish Balay           }
273a7e14dcfSSatish Balay         }
274a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
275a7e14dcfSSatish Balay       }
276a7e14dcfSSatish Balay 
277a7e14dcfSSatish Balay       if (fmin < f) {
278a7e14dcfSSatish Balay         f = fmin;
279a7e14dcfSSatish Balay         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
280a7e14dcfSSatish Balay         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
281a7e14dcfSSatish Balay 
282a9603a14SPatrick Farrell         ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
283a7e14dcfSSatish Balay 
28453506e15SBarry Smith         if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
285a7e14dcfSSatish Balay         needH = 1;
286a7e14dcfSSatish Balay 
2878931d482SJason Sarich         ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
288a7e14dcfSSatish Balay         if (reason != TAO_CONTINUE_ITERATING) {
289a7e14dcfSSatish Balay           PetscFunctionReturn(0);
290a7e14dcfSSatish Balay         }
291a7e14dcfSSatish Balay       }
292a7e14dcfSSatish Balay     }
293a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay     /*  Modify the radius if it is too large or small */
296a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tr->min_radius);
297a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
298a7e14dcfSSatish Balay     break;
299a7e14dcfSSatish Balay 
300a7e14dcfSSatish Balay   default:
301a7e14dcfSSatish Balay     /*  Norm of the first direction will initialize radius */
302a7e14dcfSSatish Balay     tao->trust = 0.0;
303a7e14dcfSSatish Balay     break;
304a7e14dcfSSatish Balay   }
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
307a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
308a7e14dcfSSatish Balay      since the function value may have decreased */
309a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type) {
310a7e14dcfSSatish Balay     if (f != 0.0) {
311a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
312a7e14dcfSSatish Balay     }
313a7e14dcfSSatish Balay     else {
314a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
315a7e14dcfSSatish Balay     }
316a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(tr->M,delta);CHKERRQ(ierr);
317a7e14dcfSSatish Balay   }
318a7e14dcfSSatish Balay 
319a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
320a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
3218931d482SJason Sarich     ++tao->niter;
322ae93cb3cSJason Sarich     tao->ksp_its=0;
323a7e14dcfSSatish Balay     /* Compute the Hessian */
324a7e14dcfSSatish Balay     if (needH) {
325ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
326a7e14dcfSSatish Balay       needH = 0;
327a7e14dcfSSatish Balay     }
328a7e14dcfSSatish Balay 
329a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type) {
330a7e14dcfSSatish Balay       if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) {
331a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner */
332a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, tr->Diag);CHKERRQ(ierr);
333a7e14dcfSSatish Balay         ierr = VecAbs(tr->Diag);CHKERRQ(ierr);
334a7e14dcfSSatish Balay         ierr = VecReciprocal(tr->Diag);CHKERRQ(ierr);
335a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(tr->M,tr->Diag);CHKERRQ(ierr);
336a7e14dcfSSatish Balay       }
337a7e14dcfSSatish Balay 
338a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
339a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
340a7e14dcfSSatish Balay       ++bfgsUpdates;
341a7e14dcfSSatish Balay     }
342a7e14dcfSSatish Balay 
343a7e14dcfSSatish Balay     while (reason == TAO_CONTINUE_ITERATING) {
34423ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
345a7e14dcfSSatish Balay 
346a7e14dcfSSatish Balay       /* Solve the trust region subproblem */
347*ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
348a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
349a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
350a7e14dcfSSatish Balay       tao->ksp_its+=its;
351ae93cb3cSJason Sarich       tao->ksp_tot_its+=its;
352*ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
353a7e14dcfSSatish Balay 
354a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
355a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
356a7e14dcfSSatish Balay         if (norm_d > 0.0) {
357a7e14dcfSSatish Balay           tao->trust = norm_d;
358a7e14dcfSSatish Balay 
359a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
360a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
361a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
362a7e14dcfSSatish Balay         }
363a7e14dcfSSatish Balay         else {
364a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
365a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
366a7e14dcfSSatish Balay           tao->trust = tao->trust0;
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
369a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
370a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
371a7e14dcfSSatish Balay 
372*ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
373a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
374a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
375a7e14dcfSSatish Balay           tao->ksp_its+=its;
3762d9aa51bSJason Sarich           tao->ksp_tot_its+=its;
377*ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
378a7e14dcfSSatish Balay 
37953506e15SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
380a7e14dcfSSatish Balay         }
381a7e14dcfSSatish Balay       }
382a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
383a7e14dcfSSatish Balay       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
384a7e14dcfSSatish Balay       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&
385a7e14dcfSSatish Balay           (NTR_PC_BFGS == tr->pc_type) && (bfgsUpdates > 1)) {
386a7e14dcfSSatish Balay         /* Preconditioner is numerically indefinite; reset the
387a7e14dcfSSatish Balay            approximate if using BFGS preconditioning. */
388a7e14dcfSSatish Balay 
389a7e14dcfSSatish Balay         if (f != 0.0) {
390a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
391a7e14dcfSSatish Balay         }
392a7e14dcfSSatish Balay         else {
393a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
394a7e14dcfSSatish Balay         }
395a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(tr->M, delta);CHKERRQ(ierr);
396a7e14dcfSSatish Balay         ierr = MatLMVMReset(tr->M);CHKERRQ(ierr);
397a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
398a7e14dcfSSatish Balay         bfgsUpdates = 1;
399a7e14dcfSSatish Balay       }
400a7e14dcfSSatish Balay 
401a7e14dcfSSatish Balay       if (NTR_UPDATE_REDUCTION == tr->update_type) {
402a7e14dcfSSatish Balay         /* Get predicted reduction */
403*ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
404a7e14dcfSSatish Balay         if (prered >= 0.0) {
405a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
406a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
407a7e14dcfSSatish Balay              be rejected! */
408a7e14dcfSSatish Balay           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
409a7e14dcfSSatish Balay         }
410a7e14dcfSSatish Balay         else {
411a7e14dcfSSatish Balay           /* Compute trial step and function value */
412a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
413a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
414a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
415a7e14dcfSSatish Balay 
416a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
417a7e14dcfSSatish Balay             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
418a7e14dcfSSatish Balay           } else {
419a7e14dcfSSatish Balay             /* Compute and actual reduction */
420a7e14dcfSSatish Balay             actred = f - ftrial;
421a7e14dcfSSatish Balay             prered = -prered;
422a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
423a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
424a7e14dcfSSatish Balay               kappa = 1.0;
425a7e14dcfSSatish Balay             }
426a7e14dcfSSatish Balay             else {
427a7e14dcfSSatish Balay               kappa = actred / prered;
428a7e14dcfSSatish Balay             }
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay             /* Accept or reject the step and update radius */
431a7e14dcfSSatish Balay             if (kappa < tr->eta1) {
432a7e14dcfSSatish Balay               /* Reject the step */
433a7e14dcfSSatish Balay               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
434a7e14dcfSSatish Balay             }
435a7e14dcfSSatish Balay             else {
436a7e14dcfSSatish Balay               /* Accept the step */
437a7e14dcfSSatish Balay               if (kappa < tr->eta2) {
438a7e14dcfSSatish Balay                 /* Marginal bad step */
439a7e14dcfSSatish Balay                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
440a7e14dcfSSatish Balay               }
441a7e14dcfSSatish Balay               else if (kappa < tr->eta3) {
442a7e14dcfSSatish Balay                 /* Reasonable step */
443a7e14dcfSSatish Balay                 tao->trust = tr->alpha3 * tao->trust;
444a7e14dcfSSatish Balay               }
445a7e14dcfSSatish Balay               else if (kappa < tr->eta4) {
446a7e14dcfSSatish Balay                 /* Good step */
447a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
448a7e14dcfSSatish Balay               }
449a7e14dcfSSatish Balay               else {
450a7e14dcfSSatish Balay                 /* Very good step */
451a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
452a7e14dcfSSatish Balay               }
453a7e14dcfSSatish Balay               break;
454a7e14dcfSSatish Balay             }
455a7e14dcfSSatish Balay           }
456a7e14dcfSSatish Balay         }
457a7e14dcfSSatish Balay       }
458a7e14dcfSSatish Balay       else {
459a7e14dcfSSatish Balay         /* Get predicted reduction */
460*ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
461a7e14dcfSSatish Balay         if (prered >= 0.0) {
462a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
463a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
464a7e14dcfSSatish Balay              be rejected! */
465a7e14dcfSSatish Balay           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
466a7e14dcfSSatish Balay         }
467a7e14dcfSSatish Balay         else {
468a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
469a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
470a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
471a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
472a7e14dcfSSatish Balay             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
473a7e14dcfSSatish Balay           }
474a7e14dcfSSatish Balay           else {
475a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
476a7e14dcfSSatish Balay             actred = f - ftrial;
477a7e14dcfSSatish Balay             prered = -prered;
478a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
479a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
480a7e14dcfSSatish Balay               kappa = 1.0;
481a7e14dcfSSatish Balay             }
482a7e14dcfSSatish Balay             else {
483a7e14dcfSSatish Balay               kappa = actred / prered;
484a7e14dcfSSatish Balay             }
485a7e14dcfSSatish Balay 
486a7e14dcfSSatish Balay             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
487a7e14dcfSSatish Balay             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
488a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
489a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
490a7e14dcfSSatish Balay 
491a7e14dcfSSatish Balay             if (kappa >= 1.0 - tr->mu1) {
492a7e14dcfSSatish Balay               /* Great agreement; accept step and update radius */
493a7e14dcfSSatish Balay               if (tau_max < 1.0) {
494a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
495a7e14dcfSSatish Balay               }
496a7e14dcfSSatish Balay               else if (tau_max > tr->gamma4) {
497a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
498a7e14dcfSSatish Balay               }
499a7e14dcfSSatish Balay               else {
500a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
501a7e14dcfSSatish Balay               }
502a7e14dcfSSatish Balay               break;
503a7e14dcfSSatish Balay             }
504a7e14dcfSSatish Balay             else if (kappa >= 1.0 - tr->mu2) {
505a7e14dcfSSatish Balay               /* Good agreement */
506a7e14dcfSSatish Balay 
507a7e14dcfSSatish Balay               if (tau_max < tr->gamma2) {
508a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
509a7e14dcfSSatish Balay               }
510a7e14dcfSSatish Balay               else if (tau_max > tr->gamma3) {
511a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
512a7e14dcfSSatish Balay               }
513a7e14dcfSSatish Balay               else if (tau_max < 1.0) {
514a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
515a7e14dcfSSatish Balay               }
516a7e14dcfSSatish Balay               else {
517a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
518a7e14dcfSSatish Balay               }
519a7e14dcfSSatish Balay               break;
520a7e14dcfSSatish Balay             }
521a7e14dcfSSatish Balay             else {
522a7e14dcfSSatish Balay               /* Not good agreement */
523a7e14dcfSSatish Balay               if (tau_min > 1.0) {
524a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
525a7e14dcfSSatish Balay               }
526a7e14dcfSSatish Balay               else if (tau_max < tr->gamma1) {
527a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
528a7e14dcfSSatish Balay               }
529a7e14dcfSSatish Balay               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
530a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
531a7e14dcfSSatish Balay               }
532a7e14dcfSSatish Balay               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
533a7e14dcfSSatish Balay                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
534a7e14dcfSSatish Balay                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
535a7e14dcfSSatish Balay               }
536a7e14dcfSSatish Balay               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
537a7e14dcfSSatish Balay                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
538a7e14dcfSSatish Balay                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
539a7e14dcfSSatish Balay               }
540a7e14dcfSSatish Balay               else {
541a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
542a7e14dcfSSatish Balay               }
543a7e14dcfSSatish Balay             }
544a7e14dcfSSatish Balay           }
545a7e14dcfSSatish Balay         }
546a7e14dcfSSatish Balay       }
547a7e14dcfSSatish Balay 
548a7e14dcfSSatish Balay       /* The step computed was not good and the radius was decreased.
549a7e14dcfSSatish Balay          Monitor the radius to terminate. */
5508931d482SJason Sarich       ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
551a7e14dcfSSatish Balay     }
552a7e14dcfSSatish Balay 
553a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
554a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay     if (reason == TAO_CONTINUE_ITERATING) {
557a7e14dcfSSatish Balay       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
558a7e14dcfSSatish Balay       f = ftrial;
559302440fdSBarry Smith       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr);
560a9603a14SPatrick Farrell       ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
56153506e15SBarry Smith       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
562a7e14dcfSSatish Balay       needH = 1;
5638931d482SJason Sarich       ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
564a7e14dcfSSatish Balay     }
565a7e14dcfSSatish Balay   }
566a7e14dcfSSatish Balay   PetscFunctionReturn(0);
567a7e14dcfSSatish Balay }
568a7e14dcfSSatish Balay 
569a7e14dcfSSatish Balay /*------------------------------------------------------------*/
570a7e14dcfSSatish Balay #undef __FUNCT__
571a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTR"
572441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao)
573a7e14dcfSSatish Balay {
574a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
575a7e14dcfSSatish Balay   PetscErrorCode ierr;
576a7e14dcfSSatish Balay 
577a7e14dcfSSatish Balay   PetscFunctionBegin;
578a7e14dcfSSatish Balay 
579a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
580a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
581a7e14dcfSSatish Balay   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
582a7e14dcfSSatish Balay 
583a7e14dcfSSatish Balay   tr->Diag = 0;
584a7e14dcfSSatish Balay   tr->M = 0;
585a7e14dcfSSatish Balay 
586a7e14dcfSSatish Balay 
587a7e14dcfSSatish Balay   PetscFunctionReturn(0);
588a7e14dcfSSatish Balay }
589a7e14dcfSSatish Balay 
590a7e14dcfSSatish Balay /*------------------------------------------------------------*/
591a7e14dcfSSatish Balay #undef __FUNCT__
592a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTR"
593441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTR(Tao tao)
594a7e14dcfSSatish Balay {
595a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
596a7e14dcfSSatish Balay   PetscErrorCode ierr;
597a7e14dcfSSatish Balay 
598a7e14dcfSSatish Balay   PetscFunctionBegin;
599a7e14dcfSSatish Balay   if (tao->setupcalled) {
600a7e14dcfSSatish Balay     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
601a7e14dcfSSatish Balay   }
602a7e14dcfSSatish Balay   ierr = MatDestroy(&tr->M);CHKERRQ(ierr);
603a7e14dcfSSatish Balay   ierr = VecDestroy(&tr->Diag);CHKERRQ(ierr);
604a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
605a7e14dcfSSatish Balay   PetscFunctionReturn(0);
606a7e14dcfSSatish Balay }
607a7e14dcfSSatish Balay 
608a7e14dcfSSatish Balay /*------------------------------------------------------------*/
609a7e14dcfSSatish Balay #undef __FUNCT__
610a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTR"
6114416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(PetscOptionItems *PetscOptionsObject,Tao tao)
612a7e14dcfSSatish Balay {
613a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
614a7e14dcfSSatish Balay   PetscErrorCode ierr;
615a7e14dcfSSatish Balay 
616a7e14dcfSSatish Balay   PetscFunctionBegin;
6171a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
61894ae4db5SBarry 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);
61994ae4db5SBarry 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);
62094ae4db5SBarry 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);
62194ae4db5SBarry 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);
62294ae4db5SBarry 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);
62394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr);
62494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr);
62594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr);
62694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr);
62794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr);
62894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr);
62994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr);
63094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr);
63194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr);
63294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr);
63394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr);
63494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr);
63594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr);
63694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr);
63794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr);
63894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr);
63994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr);
64094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr);
64194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr);
64294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr);
64394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr);
64494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr);
64594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr);
64694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr);
64794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr);
64894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr);
649a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
650a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
651a7e14dcfSSatish Balay   PetscFunctionReturn(0);
652a7e14dcfSSatish Balay }
653a7e14dcfSSatish Balay 
654a7e14dcfSSatish Balay /*------------------------------------------------------------*/
655a7e14dcfSSatish Balay #undef __FUNCT__
656a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTR"
657441846f8SBarry Smith static PetscErrorCode TaoView_NTR(Tao tao, PetscViewer viewer)
658a7e14dcfSSatish Balay {
659a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
660a7e14dcfSSatish Balay   PetscErrorCode ierr;
661a7e14dcfSSatish Balay   PetscInt       nrejects;
662a7e14dcfSSatish Balay   PetscBool      isascii;
66353506e15SBarry Smith 
664a7e14dcfSSatish Balay   PetscFunctionBegin;
665a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
666a7e14dcfSSatish Balay   if (isascii) {
667a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
668a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type && tr->M) {
669a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(tr->M, &nrejects);CHKERRQ(ierr);
670a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr);
671a7e14dcfSSatish Balay     }
672a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
673a7e14dcfSSatish Balay   }
674a7e14dcfSSatish Balay   PetscFunctionReturn(0);
675a7e14dcfSSatish Balay }
676a7e14dcfSSatish Balay 
677a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6781522df2eSJason Sarich /*MC
6791522df2eSJason Sarich   TAONTR - Newton's method with trust region for unconstrained minimization.
6801522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system.
6811522df2eSJason Sarich   NTR expects a KSP solver with a trust region radius.
6821522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
6831522df2eSJason Sarich 
6841522df2eSJason Sarich   Options Database Keys:
6851522df2eSJason Sarich + -tao_ntr_ksp_type - "nash","stcg","gltr"
6861522df2eSJason Sarich . -tao_ntr_pc_type - "none","ahess","bfgs","petsc"
6871522df2eSJason Sarich . -tao_ntr_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs"
6881522df2eSJason Sarich . -tao_ntr_init_type - "constant","direction","interpolation"
6891522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation"
6901522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius
6911522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius
6921522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
6931522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor
6941522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor
6951522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor
6961522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor
6971522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor
6981522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor
6991522df2eSJason Sarich . -tao_ntr_theta_i - thetha1 interpolation init factor
7001522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor
7011522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor
7021522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor
7031522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor
7041522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor
7051522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor
7061522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor
7071522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
7081522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
7091522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update
7101522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update
7111522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update
7121522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update
7131522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update
7141522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update
7151522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update
7161522df2eSJason Sarich 
7171eb8069cSJason Sarich   Level: beginner
7181522df2eSJason Sarich M*/
7191522df2eSJason Sarich 
720a7e14dcfSSatish Balay #undef __FUNCT__
721a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTR"
722728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
723a7e14dcfSSatish Balay {
724a7e14dcfSSatish Balay   TAO_NTR *tr;
725a7e14dcfSSatish Balay   PetscErrorCode ierr;
726a7e14dcfSSatish Balay 
727a7e14dcfSSatish Balay   PetscFunctionBegin;
728a7e14dcfSSatish Balay 
7293c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
730a7e14dcfSSatish Balay 
731a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTR;
732a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTR;
733a7e14dcfSSatish Balay   tao->ops->view = TaoView_NTR;
734a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
735a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTR;
736a7e14dcfSSatish Balay 
7376552cf8aSJason Sarich   /* Override default settings (unless already changed) */
7386552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
7396552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
740a7e14dcfSSatish Balay   tao->data = (void*)tr;
741a7e14dcfSSatish Balay 
742a7e14dcfSSatish Balay   /*  Standard trust region update parameters */
743a7e14dcfSSatish Balay   tr->eta1 = 1.0e-4;
744a7e14dcfSSatish Balay   tr->eta2 = 0.25;
745a7e14dcfSSatish Balay   tr->eta3 = 0.50;
746a7e14dcfSSatish Balay   tr->eta4 = 0.90;
747a7e14dcfSSatish Balay 
748a7e14dcfSSatish Balay   tr->alpha1 = 0.25;
749a7e14dcfSSatish Balay   tr->alpha2 = 0.50;
750a7e14dcfSSatish Balay   tr->alpha3 = 1.00;
751a7e14dcfSSatish Balay   tr->alpha4 = 2.00;
752a7e14dcfSSatish Balay   tr->alpha5 = 4.00;
753a7e14dcfSSatish Balay 
754a7e14dcfSSatish Balay   /*  Interpolation parameters */
755a7e14dcfSSatish Balay   tr->mu1_i = 0.35;
756a7e14dcfSSatish Balay   tr->mu2_i = 0.50;
757a7e14dcfSSatish Balay 
758a7e14dcfSSatish Balay   tr->gamma1_i = 0.0625;
759a7e14dcfSSatish Balay   tr->gamma2_i = 0.50;
760a7e14dcfSSatish Balay   tr->gamma3_i = 2.00;
761a7e14dcfSSatish Balay   tr->gamma4_i = 5.00;
762a7e14dcfSSatish Balay 
763a7e14dcfSSatish Balay   tr->theta_i = 0.25;
764a7e14dcfSSatish Balay 
765a7e14dcfSSatish Balay   /*  Interpolation trust region update parameters */
766a7e14dcfSSatish Balay   tr->mu1 = 0.10;
767a7e14dcfSSatish Balay   tr->mu2 = 0.50;
768a7e14dcfSSatish Balay 
769a7e14dcfSSatish Balay   tr->gamma1 = 0.25;
770a7e14dcfSSatish Balay   tr->gamma2 = 0.50;
771a7e14dcfSSatish Balay   tr->gamma3 = 2.00;
772a7e14dcfSSatish Balay   tr->gamma4 = 4.00;
773a7e14dcfSSatish Balay 
774a7e14dcfSSatish Balay   tr->theta = 0.05;
775a7e14dcfSSatish Balay 
776a7e14dcfSSatish Balay   tr->min_radius = 1.0e-10;
777a7e14dcfSSatish Balay   tr->max_radius = 1.0e10;
778a7e14dcfSSatish Balay   tr->epsilon = 1.0e-6;
779a7e14dcfSSatish Balay 
780a7e14dcfSSatish Balay   tr->ksp_type        = NTR_KSP_STCG;
781a7e14dcfSSatish Balay   tr->pc_type         = NTR_PC_BFGS;
782a7e14dcfSSatish Balay   tr->bfgs_scale_type = BFGS_SCALE_AHESS;
783a7e14dcfSSatish Balay   tr->init_type       = NTR_INIT_INTERPOLATION;
784a7e14dcfSSatish Balay   tr->update_type     = NTR_UPDATE_REDUCTION;
785a7e14dcfSSatish Balay 
786a7e14dcfSSatish Balay 
787a7e14dcfSSatish Balay   /* Set linear solver to default for trust region */
788a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
7895d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
790a7e14dcfSSatish Balay   PetscFunctionReturn(0);
791a7e14dcfSSatish Balay }
792a7e14dcfSSatish Balay 
793a7e14dcfSSatish Balay 
794a7e14dcfSSatish Balay #undef __FUNCT__
795a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
796a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
797a7e14dcfSSatish Balay {
798a7e14dcfSSatish Balay     PetscErrorCode ierr;
799a7e14dcfSSatish Balay     Mat M;
800a7e14dcfSSatish Balay     PetscFunctionBegin;
801a7e14dcfSSatish Balay     PetscValidHeaderSpecific(pc,PC_CLASSID,1);
802a7e14dcfSSatish Balay     PetscValidHeaderSpecific(b,VEC_CLASSID,2);
803a7e14dcfSSatish Balay     PetscValidHeaderSpecific(x,VEC_CLASSID,3);
804a7e14dcfSSatish Balay     ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
805a7e14dcfSSatish Balay     ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
806a7e14dcfSSatish Balay     PetscFunctionReturn(0);
807a7e14dcfSSatish Balay }
808