xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision 3c859ba3a04a72e8efcb87bc7ffc046a6cbab413)
1fb90e4d1STodd Munson #include <../src/tao/unconstrained/impls/ntr/ntrimpl.h>
2a7e14dcfSSatish Balay 
3aaa7dc30SBarry Smith #include <petscksp.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay #define NTR_INIT_CONSTANT         0
6a7e14dcfSSatish Balay #define NTR_INIT_DIRECTION        1
7a7e14dcfSSatish Balay #define NTR_INIT_INTERPOLATION    2
8a7e14dcfSSatish Balay #define NTR_INIT_TYPES            3
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay #define NTR_UPDATE_REDUCTION      0
11a7e14dcfSSatish Balay #define NTR_UPDATE_INTERPOLATION  1
12a7e14dcfSSatish Balay #define NTR_UPDATE_TYPES          2
13a7e14dcfSSatish Balay 
1453506e15SBarry Smith static const char *NTR_INIT[64] = {"constant","direction","interpolation"};
15a7e14dcfSSatish Balay 
1653506e15SBarry Smith static const char *NTR_UPDATE[64] = {"reduction","interpolation"};
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay /*
19a7e14dcfSSatish Balay    TaoSolve_NTR - Implements Newton's Method with a trust region approach
20a7e14dcfSSatish Balay    for solving unconstrained minimization problems.
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay    The basic algorithm is taken from MINPACK-2 (dstrn).
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay    TaoSolve_NTR computes a local minimizer of a twice differentiable function
25a7e14dcfSSatish Balay    f by applying a trust region variant of Newton's method.  At each stage
26a7e14dcfSSatish Balay    of the algorithm, we use the prconditioned conjugate gradient method to
27a7e14dcfSSatish Balay    determine an approximate minimizer of the quadratic equation
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay         q(s) = <s, Hs + g>
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay    subject to the trust region constraint
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay         || s ||_M <= radius,
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay    where radius is the trust region radius and M is a symmetric positive
36a7e14dcfSSatish Balay    definite matrix (the preconditioner).  Here g is the gradient and H
37a7e14dcfSSatish Balay    is the Hessian matrix.
38a7e14dcfSSatish Balay 
3905de396fSBarry Smith    Note:  TaoSolve_NTR MUST use the iterative solver KSPNASH, KSPSTCG,
4005de396fSBarry Smith           or KSPGLTR.  Thus, we set KSPNASH, KSPSTCG, or KSPGLTR in this
41a7e14dcfSSatish Balay           routine regardless of what the user may have previously specified.
42a7e14dcfSSatish Balay */
43441846f8SBarry Smith static PetscErrorCode TaoSolve_NTR(Tao tao)
44a7e14dcfSSatish Balay {
45a7e14dcfSSatish Balay   TAO_NTR            *tr = (TAO_NTR *)tao->data;
46fb90e4d1STodd Munson   KSPType            ksp_type;
470ad3a497SAlp Dener   PetscBool          is_nash,is_stcg,is_gltr,is_bfgs,is_jacobi,is_symmetric,sym_set;
48a7e14dcfSSatish Balay   KSPConvergedReason ksp_reason;
49fb90e4d1STodd Munson   PC                 pc;
50a7e14dcfSSatish Balay   PetscReal          fmin, ftrial, prered, actred, kappa, sigma, beta;
51a7e14dcfSSatish Balay   PetscReal          tau, tau_1, tau_2, tau_max, tau_min, max_radius;
52a7e14dcfSSatish Balay   PetscReal          f, gnorm;
53a7e14dcfSSatish Balay 
54a7e14dcfSSatish Balay   PetscReal          norm_d;
55a7e14dcfSSatish Balay   PetscErrorCode     ierr;
56a7e14dcfSSatish Balay   PetscInt           bfgsUpdates = 0;
57a7e14dcfSSatish Balay   PetscInt           needH;
58a7e14dcfSSatish Balay 
59a7e14dcfSSatish Balay   PetscInt           i_max = 5;
60a7e14dcfSSatish Balay   PetscInt           j_max = 1;
61a7e14dcfSSatish Balay   PetscInt           i, j, N, n, its;
62a7e14dcfSSatish Balay 
63a7e14dcfSSatish Balay   PetscFunctionBegin;
64a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
659dcef436SStefano Zampini     ierr = PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr);
66a7e14dcfSSatish Balay   }
67a7e14dcfSSatish Balay 
68fb90e4d1STodd Munson   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
6905de396fSBarry Smith   ierr = PetscStrcmp(ksp_type,KSPNASH,&is_nash);CHKERRQ(ierr);
7005de396fSBarry Smith   ierr = PetscStrcmp(ksp_type,KSPSTCG,&is_stcg);CHKERRQ(ierr);
7105de396fSBarry Smith   ierr = PetscStrcmp(ksp_type,KSPGLTR,&is_gltr);CHKERRQ(ierr);
72*3c859ba3SBarry Smith   PetscCheck(is_nash || is_stcg || is_gltr,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"TAO_NTR requires nash, stcg, or gltr for the KSP");
73a7e14dcfSSatish Balay 
74fb90e4d1STodd Munson   /* Initialize the radius and modify if it is too large or small */
75fb90e4d1STodd Munson   tao->trust = tao->trust0;
76a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tr->min_radius);
77a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tr->max_radius);
78a7e14dcfSSatish Balay 
790c51296cSAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
800c51296cSAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
810c51296cSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr);
820c51296cSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr);
830c51296cSAlp Dener   if (is_bfgs) {
840c51296cSAlp Dener     tr->bfgs_pre = pc;
850c51296cSAlp Dener     ierr = PCLMVMGetMatLMVM(tr->bfgs_pre, &tr->M);CHKERRQ(ierr);
86a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
87a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
880c51296cSAlp Dener     ierr = MatSetSizes(tr->M, n, n, N, N);CHKERRQ(ierr);
89cd929ea3SAlp Dener     ierr = MatLMVMAllocate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
900ad3a497SAlp Dener     ierr = MatIsSymmetricKnown(tr->M, &sym_set, &is_symmetric);CHKERRQ(ierr);
91*3c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
920c51296cSAlp Dener   } else if (is_jacobi) {
930c51296cSAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
94a7e14dcfSSatish Balay   }
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay   /* Check convergence criteria */
97a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
98a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
99*3c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER,"User provided compute function generated Inf or NaN");
100a7e14dcfSSatish Balay   needH = 1;
101a7e14dcfSSatish Balay 
1023ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1033ecd9318SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1043ecd9318SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,1.0);CHKERRQ(ierr);
1053ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1063ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
107a7e14dcfSSatish Balay 
108a7e14dcfSSatish Balay   /*  Initialize trust-region radius */
109a7e14dcfSSatish Balay   switch (tr->init_type) {
110a7e14dcfSSatish Balay   case NTR_INIT_CONSTANT:
111a7e14dcfSSatish Balay     /*  Use the initial radius specified */
112a7e14dcfSSatish Balay     break;
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay   case NTR_INIT_INTERPOLATION:
115a7e14dcfSSatish Balay     /*  Use the initial radius specified */
116a7e14dcfSSatish Balay     max_radius = 0.0;
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
119a7e14dcfSSatish Balay       fmin = f;
120a7e14dcfSSatish Balay       sigma = 0.0;
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay       if (needH) {
123ffad9901SBarry Smith         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
124a7e14dcfSSatish Balay         needH = 0;
125a7e14dcfSSatish Balay       }
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
130a7e14dcfSSatish Balay         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
131a7e14dcfSSatish Balay         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
134a7e14dcfSSatish Balay           tau = tr->gamma1_i;
135a7e14dcfSSatish Balay         }
136a7e14dcfSSatish Balay         else {
137a7e14dcfSSatish Balay           if (ftrial < fmin) {
138a7e14dcfSSatish Balay             fmin = ftrial;
139a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
140a7e14dcfSSatish Balay           }
141a7e14dcfSSatish Balay 
142a7e14dcfSSatish Balay           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
143a7e14dcfSSatish Balay           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
146a7e14dcfSSatish Balay           actred = f - ftrial;
147a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
148a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tr->epsilon)) {
149a7e14dcfSSatish Balay             kappa = 1.0;
150a7e14dcfSSatish Balay           }
151a7e14dcfSSatish Balay           else {
152a7e14dcfSSatish Balay             kappa = actred / prered;
153a7e14dcfSSatish Balay           }
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
156a7e14dcfSSatish Balay           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
157a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
158a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
159a7e14dcfSSatish Balay 
16018cfbf8eSSatish Balay           if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tr->mu1_i) {
161a7e14dcfSSatish Balay             /*  Great agreement */
162a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay             if (tau_max < 1.0) {
165a7e14dcfSSatish Balay               tau = tr->gamma3_i;
166a7e14dcfSSatish Balay             }
167a7e14dcfSSatish Balay             else if (tau_max > tr->gamma4_i) {
168a7e14dcfSSatish Balay               tau = tr->gamma4_i;
169a7e14dcfSSatish Balay             }
170a7e14dcfSSatish Balay             else {
171a7e14dcfSSatish Balay               tau = tau_max;
172a7e14dcfSSatish Balay             }
173a7e14dcfSSatish Balay           }
17418cfbf8eSSatish Balay           else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tr->mu2_i) {
175a7e14dcfSSatish Balay             /*  Good agreement */
176a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay             if (tau_max < tr->gamma2_i) {
179a7e14dcfSSatish Balay               tau = tr->gamma2_i;
180a7e14dcfSSatish Balay             }
181a7e14dcfSSatish Balay             else if (tau_max > tr->gamma3_i) {
182a7e14dcfSSatish Balay               tau = tr->gamma3_i;
183a7e14dcfSSatish Balay             }
184a7e14dcfSSatish Balay             else {
185a7e14dcfSSatish Balay               tau = tau_max;
186a7e14dcfSSatish Balay             }
187a7e14dcfSSatish Balay           }
188a7e14dcfSSatish Balay           else {
189a7e14dcfSSatish Balay             /*  Not good agreement */
190a7e14dcfSSatish Balay             if (tau_min > 1.0) {
191a7e14dcfSSatish Balay               tau = tr->gamma2_i;
192a7e14dcfSSatish Balay             }
193a7e14dcfSSatish Balay             else if (tau_max < tr->gamma1_i) {
194a7e14dcfSSatish Balay               tau = tr->gamma1_i;
195a7e14dcfSSatish Balay             }
196a7e14dcfSSatish Balay             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
197a7e14dcfSSatish Balay               tau = tr->gamma1_i;
198a7e14dcfSSatish Balay             }
199a7e14dcfSSatish Balay             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
200a7e14dcfSSatish Balay                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
201a7e14dcfSSatish Balay               tau = tau_1;
202a7e14dcfSSatish Balay             }
203a7e14dcfSSatish Balay             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
204a7e14dcfSSatish Balay                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
205a7e14dcfSSatish Balay               tau = tau_2;
206a7e14dcfSSatish Balay             }
207a7e14dcfSSatish Balay             else {
208a7e14dcfSSatish Balay               tau = tau_max;
209a7e14dcfSSatish Balay             }
210a7e14dcfSSatish Balay           }
211a7e14dcfSSatish Balay         }
212a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
213a7e14dcfSSatish Balay       }
214a7e14dcfSSatish Balay 
215a7e14dcfSSatish Balay       if (fmin < f) {
216a7e14dcfSSatish Balay         f = fmin;
217a7e14dcfSSatish Balay         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
218a7e14dcfSSatish Balay         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
219a7e14dcfSSatish Balay 
220a9603a14SPatrick Farrell         ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
221*3c859ba3SBarry Smith         PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
222a7e14dcfSSatish Balay         needH = 1;
223a7e14dcfSSatish Balay 
2243ecd9318SAlp Dener         ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
2253ecd9318SAlp Dener         ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,1.0);CHKERRQ(ierr);
2263ecd9318SAlp Dener         ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2273ecd9318SAlp Dener         if (tao->reason != TAO_CONTINUE_ITERATING) {
228a7e14dcfSSatish Balay           PetscFunctionReturn(0);
229a7e14dcfSSatish Balay         }
230a7e14dcfSSatish Balay       }
231a7e14dcfSSatish Balay     }
232a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
233a7e14dcfSSatish Balay 
234a7e14dcfSSatish Balay     /*  Modify the radius if it is too large or small */
235a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tr->min_radius);
236a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
237a7e14dcfSSatish Balay     break;
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay   default:
240a7e14dcfSSatish Balay     /*  Norm of the first direction will initialize radius */
241a7e14dcfSSatish Balay     tao->trust = 0.0;
242a7e14dcfSSatish Balay     break;
243a7e14dcfSSatish Balay   }
244a7e14dcfSSatish Balay 
245a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2463ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
247e1e80dc8SAlp Dener     /* Call general purpose update function */
248e1e80dc8SAlp Dener     if (tao->ops->update) {
2498fcddce6SStefano Zampini       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
250e1e80dc8SAlp Dener     }
2518931d482SJason Sarich     ++tao->niter;
252ae93cb3cSJason Sarich     tao->ksp_its=0;
253a7e14dcfSSatish Balay     /* Compute the Hessian */
254a7e14dcfSSatish Balay     if (needH) {
255ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
256a7e14dcfSSatish Balay       needH = 0;
257a7e14dcfSSatish Balay     }
258a7e14dcfSSatish Balay 
2590c51296cSAlp Dener     if (tr->bfgs_pre) {
260a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
261a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
262a7e14dcfSSatish Balay       ++bfgsUpdates;
263a7e14dcfSSatish Balay     }
264a7e14dcfSSatish Balay 
2653ecd9318SAlp Dener     while (tao->reason == TAO_CONTINUE_ITERATING) {
26623ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay       /* Solve the trust region subproblem */
269ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
270a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
271a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
272a7e14dcfSSatish Balay       tao->ksp_its+=its;
273ae93cb3cSJason Sarich       tao->ksp_tot_its+=its;
274ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
275a7e14dcfSSatish Balay 
276a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
277a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
278a7e14dcfSSatish Balay         if (norm_d > 0.0) {
279a7e14dcfSSatish Balay           tao->trust = norm_d;
280a7e14dcfSSatish Balay 
281a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
282a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
283a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
284a7e14dcfSSatish Balay         }
285a7e14dcfSSatish Balay         else {
286a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
287a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
288a7e14dcfSSatish Balay           tao->trust = tao->trust0;
289a7e14dcfSSatish Balay 
290a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
291a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
292a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
293a7e14dcfSSatish Balay 
294ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
295a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
296a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
297a7e14dcfSSatish Balay           tao->ksp_its+=its;
2982d9aa51bSJason Sarich           tao->ksp_tot_its+=its;
299ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
300a7e14dcfSSatish Balay 
301*3c859ba3SBarry Smith           PetscCheck(norm_d != 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero");
302a7e14dcfSSatish Balay         }
303a7e14dcfSSatish Balay       }
304a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
305a7e14dcfSSatish Balay       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
3060c51296cSAlp Dener       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (tr->bfgs_pre)) {
307a7e14dcfSSatish Balay         /* Preconditioner is numerically indefinite; reset the
308a7e14dcfSSatish Balay            approximate if using BFGS preconditioning. */
309cd929ea3SAlp Dener         ierr = MatLMVMReset(tr->M, PETSC_FALSE);CHKERRQ(ierr);
310a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
311a7e14dcfSSatish Balay         bfgsUpdates = 1;
312a7e14dcfSSatish Balay       }
313a7e14dcfSSatish Balay 
314a7e14dcfSSatish Balay       if (NTR_UPDATE_REDUCTION == tr->update_type) {
315a7e14dcfSSatish Balay         /* Get predicted reduction */
316ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
317a7e14dcfSSatish Balay         if (prered >= 0.0) {
318a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
319a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
320a7e14dcfSSatish Balay              be rejected! */
321a7e14dcfSSatish Balay           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
322a7e14dcfSSatish Balay         }
323a7e14dcfSSatish Balay         else {
324a7e14dcfSSatish Balay           /* Compute trial step and function value */
325a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
326a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
327a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
328a7e14dcfSSatish Balay 
329a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
330a7e14dcfSSatish Balay             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
331a7e14dcfSSatish Balay           } else {
332a7e14dcfSSatish Balay             /* Compute and actual reduction */
333a7e14dcfSSatish Balay             actred = f - ftrial;
334a7e14dcfSSatish Balay             prered = -prered;
335a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
336a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
337a7e14dcfSSatish Balay               kappa = 1.0;
338a7e14dcfSSatish Balay             }
339a7e14dcfSSatish Balay             else {
340a7e14dcfSSatish Balay               kappa = actred / prered;
341a7e14dcfSSatish Balay             }
342a7e14dcfSSatish Balay 
343a7e14dcfSSatish Balay             /* Accept or reject the step and update radius */
344a7e14dcfSSatish Balay             if (kappa < tr->eta1) {
345a7e14dcfSSatish Balay               /* Reject the step */
346a7e14dcfSSatish Balay               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
347a7e14dcfSSatish Balay             }
348a7e14dcfSSatish Balay             else {
349a7e14dcfSSatish Balay               /* Accept the step */
350a7e14dcfSSatish Balay               if (kappa < tr->eta2) {
351a7e14dcfSSatish Balay                 /* Marginal bad step */
352a7e14dcfSSatish Balay                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
353a7e14dcfSSatish Balay               }
354a7e14dcfSSatish Balay               else if (kappa < tr->eta3) {
355a7e14dcfSSatish Balay                 /* Reasonable step */
356a7e14dcfSSatish Balay                 tao->trust = tr->alpha3 * tao->trust;
357a7e14dcfSSatish Balay               }
358a7e14dcfSSatish Balay               else if (kappa < tr->eta4) {
359a7e14dcfSSatish Balay                 /* Good step */
360a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
361a7e14dcfSSatish Balay               }
362a7e14dcfSSatish Balay               else {
363a7e14dcfSSatish Balay                 /* Very good step */
364a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
365a7e14dcfSSatish Balay               }
366a7e14dcfSSatish Balay               break;
367a7e14dcfSSatish Balay             }
368a7e14dcfSSatish Balay           }
369a7e14dcfSSatish Balay         }
370a7e14dcfSSatish Balay       }
371a7e14dcfSSatish Balay       else {
372a7e14dcfSSatish Balay         /* Get predicted reduction */
373ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
374a7e14dcfSSatish Balay         if (prered >= 0.0) {
375a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
376a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
377a7e14dcfSSatish Balay              be rejected! */
378a7e14dcfSSatish Balay           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
379a7e14dcfSSatish Balay         }
380a7e14dcfSSatish Balay         else {
381a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
382a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
383a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
384a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
385a7e14dcfSSatish Balay             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
386a7e14dcfSSatish Balay           }
387a7e14dcfSSatish Balay           else {
388a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
389a7e14dcfSSatish Balay             actred = f - ftrial;
390a7e14dcfSSatish Balay             prered = -prered;
391a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
392a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
393a7e14dcfSSatish Balay               kappa = 1.0;
394a7e14dcfSSatish Balay             }
395a7e14dcfSSatish Balay             else {
396a7e14dcfSSatish Balay               kappa = actred / prered;
397a7e14dcfSSatish Balay             }
398a7e14dcfSSatish Balay 
399a7e14dcfSSatish Balay             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
400a7e14dcfSSatish Balay             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
401a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
402a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
403a7e14dcfSSatish Balay 
404a7e14dcfSSatish Balay             if (kappa >= 1.0 - tr->mu1) {
405a7e14dcfSSatish Balay               /* Great agreement; accept step and update radius */
406a7e14dcfSSatish Balay               if (tau_max < 1.0) {
407a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
408a7e14dcfSSatish Balay               }
409a7e14dcfSSatish Balay               else if (tau_max > tr->gamma4) {
410a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
411a7e14dcfSSatish Balay               }
412a7e14dcfSSatish Balay               else {
413a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
414a7e14dcfSSatish Balay               }
415a7e14dcfSSatish Balay               break;
416a7e14dcfSSatish Balay             }
417a7e14dcfSSatish Balay             else if (kappa >= 1.0 - tr->mu2) {
418a7e14dcfSSatish Balay               /* Good agreement */
419a7e14dcfSSatish Balay 
420a7e14dcfSSatish Balay               if (tau_max < tr->gamma2) {
421a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
422a7e14dcfSSatish Balay               }
423a7e14dcfSSatish Balay               else if (tau_max > tr->gamma3) {
424a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
425a7e14dcfSSatish Balay               }
426a7e14dcfSSatish Balay               else if (tau_max < 1.0) {
427a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
428a7e14dcfSSatish Balay               }
429a7e14dcfSSatish Balay               else {
430a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
431a7e14dcfSSatish Balay               }
432a7e14dcfSSatish Balay               break;
433a7e14dcfSSatish Balay             }
434a7e14dcfSSatish Balay             else {
435a7e14dcfSSatish Balay               /* Not good agreement */
436a7e14dcfSSatish Balay               if (tau_min > 1.0) {
437a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
438a7e14dcfSSatish Balay               }
439a7e14dcfSSatish Balay               else if (tau_max < tr->gamma1) {
440a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
441a7e14dcfSSatish Balay               }
442a7e14dcfSSatish Balay               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
443a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
444a7e14dcfSSatish Balay               }
445a7e14dcfSSatish Balay               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
446a7e14dcfSSatish Balay                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
447a7e14dcfSSatish Balay                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
448a7e14dcfSSatish Balay               }
449a7e14dcfSSatish Balay               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
450a7e14dcfSSatish Balay                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
451a7e14dcfSSatish Balay                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
452a7e14dcfSSatish Balay               }
453a7e14dcfSSatish Balay               else {
454a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
455a7e14dcfSSatish Balay               }
456a7e14dcfSSatish Balay             }
457a7e14dcfSSatish Balay           }
458a7e14dcfSSatish Balay         }
459a7e14dcfSSatish Balay       }
460a7e14dcfSSatish Balay 
461a7e14dcfSSatish Balay       /* The step computed was not good and the radius was decreased.
462a7e14dcfSSatish Balay          Monitor the radius to terminate. */
4633ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
4643ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,tao->trust);CHKERRQ(ierr);
4653ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
466a7e14dcfSSatish Balay     }
467a7e14dcfSSatish Balay 
468a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
469a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
470a7e14dcfSSatish Balay 
4713ecd9318SAlp Dener     if (tao->reason == TAO_CONTINUE_ITERATING) {
472a7e14dcfSSatish Balay       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
473a7e14dcfSSatish Balay       f = ftrial;
474302440fdSBarry Smith       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr);
475a9603a14SPatrick Farrell       ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
476*3c859ba3SBarry Smith       PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
477a7e14dcfSSatish Balay       needH = 1;
4783ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
4793ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,tao->trust);CHKERRQ(ierr);
4803ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
481a7e14dcfSSatish Balay     }
482a7e14dcfSSatish Balay   }
483a7e14dcfSSatish Balay   PetscFunctionReturn(0);
484a7e14dcfSSatish Balay }
485a7e14dcfSSatish Balay 
486a7e14dcfSSatish Balay /*------------------------------------------------------------*/
487441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao)
488a7e14dcfSSatish Balay {
489a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
490a7e14dcfSSatish Balay   PetscErrorCode ierr;
491a7e14dcfSSatish Balay 
492a7e14dcfSSatish Balay   PetscFunctionBegin;
493a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
494a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
495a7e14dcfSSatish Balay   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
496a7e14dcfSSatish Balay 
49783c8fe1dSLisandro Dalcin   tr->bfgs_pre = NULL;
49883c8fe1dSLisandro Dalcin   tr->M = NULL;
499a7e14dcfSSatish Balay   PetscFunctionReturn(0);
500a7e14dcfSSatish Balay }
501a7e14dcfSSatish Balay 
502a7e14dcfSSatish Balay /*------------------------------------------------------------*/
503441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTR(Tao tao)
504a7e14dcfSSatish Balay {
505a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
506a7e14dcfSSatish Balay   PetscErrorCode ierr;
507a7e14dcfSSatish Balay 
508a7e14dcfSSatish Balay   PetscFunctionBegin;
509a7e14dcfSSatish Balay   if (tao->setupcalled) {
510a7e14dcfSSatish Balay     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
511a7e14dcfSSatish Balay   }
512a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
513a7e14dcfSSatish Balay   PetscFunctionReturn(0);
514a7e14dcfSSatish Balay }
515a7e14dcfSSatish Balay 
516a7e14dcfSSatish Balay /*------------------------------------------------------------*/
5174416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(PetscOptionItems *PetscOptionsObject,Tao tao)
518a7e14dcfSSatish Balay {
519a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
520a7e14dcfSSatish Balay   PetscErrorCode ierr;
521a7e14dcfSSatish Balay 
522a7e14dcfSSatish Balay   PetscFunctionBegin;
5231a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
52494ae4db5SBarry 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);
52594ae4db5SBarry 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);
52694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr);
52794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr);
52894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr);
52994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr);
53094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr);
53194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr);
53294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr);
53394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr);
53494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr);
53594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr);
53694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr);
53794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr);
53894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr);
53994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr);
54094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr);
54194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr);
54294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr);
54394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr);
54494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr);
54594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr);
54694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr);
54794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr);
54894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr);
54994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr);
55094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr);
55194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr);
552a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
553a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
554a7e14dcfSSatish Balay   PetscFunctionReturn(0);
555a7e14dcfSSatish Balay }
556a7e14dcfSSatish Balay 
557a7e14dcfSSatish Balay /*------------------------------------------------------------*/
5581522df2eSJason Sarich /*MC
5591522df2eSJason Sarich   TAONTR - Newton's method with trust region for unconstrained minimization.
5601522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system.
5611522df2eSJason Sarich   NTR expects a KSP solver with a trust region radius.
5621522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
5631522df2eSJason Sarich 
5641522df2eSJason Sarich   Options Database Keys:
5659d0a60b2SAlp Dener + -tao_ntr_init_type - "constant","direction","interpolation"
5661522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation"
5671522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius
5681522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius
5691522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
5701522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor
5711522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor
5721522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor
5731522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor
5741522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor
5751522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor
5768966356dSPierre Jolivet . -tao_ntr_theta_i - theta1 interpolation init factor
5771522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor
5781522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor
5791522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor
5801522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor
5811522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor
5821522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor
5831522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor
5841522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
5851522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
5861522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update
5871522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update
5881522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update
5891522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update
5901522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update
5911522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update
5921522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update
5931522df2eSJason Sarich 
5941eb8069cSJason Sarich   Level: beginner
5951522df2eSJason Sarich M*/
5961522df2eSJason Sarich 
597728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
598a7e14dcfSSatish Balay {
599a7e14dcfSSatish Balay   TAO_NTR *tr;
600a7e14dcfSSatish Balay   PetscErrorCode ierr;
601a7e14dcfSSatish Balay 
602a7e14dcfSSatish Balay   PetscFunctionBegin;
603a7e14dcfSSatish Balay 
6043c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
605a7e14dcfSSatish Balay 
606a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTR;
607a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTR;
608a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
609a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTR;
610a7e14dcfSSatish Balay 
6116552cf8aSJason Sarich   /* Override default settings (unless already changed) */
6126552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
6136552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
614a7e14dcfSSatish Balay   tao->data = (void*)tr;
615a7e14dcfSSatish Balay 
616a7e14dcfSSatish Balay   /*  Standard trust region update parameters */
617a7e14dcfSSatish Balay   tr->eta1 = 1.0e-4;
618a7e14dcfSSatish Balay   tr->eta2 = 0.25;
619a7e14dcfSSatish Balay   tr->eta3 = 0.50;
620a7e14dcfSSatish Balay   tr->eta4 = 0.90;
621a7e14dcfSSatish Balay 
622a7e14dcfSSatish Balay   tr->alpha1 = 0.25;
623a7e14dcfSSatish Balay   tr->alpha2 = 0.50;
624a7e14dcfSSatish Balay   tr->alpha3 = 1.00;
625a7e14dcfSSatish Balay   tr->alpha4 = 2.00;
626a7e14dcfSSatish Balay   tr->alpha5 = 4.00;
627a7e14dcfSSatish Balay 
628a7e14dcfSSatish Balay   /*  Interpolation trust region update parameters */
629a7e14dcfSSatish Balay   tr->mu1 = 0.10;
630a7e14dcfSSatish Balay   tr->mu2 = 0.50;
631a7e14dcfSSatish Balay 
632a7e14dcfSSatish Balay   tr->gamma1 = 0.25;
633a7e14dcfSSatish Balay   tr->gamma2 = 0.50;
634a7e14dcfSSatish Balay   tr->gamma3 = 2.00;
635a7e14dcfSSatish Balay   tr->gamma4 = 4.00;
636a7e14dcfSSatish Balay 
637a7e14dcfSSatish Balay   tr->theta = 0.05;
638a7e14dcfSSatish Balay 
639fb90e4d1STodd Munson   /*  Interpolation parameters for initialization */
640fb90e4d1STodd Munson   tr->mu1_i = 0.35;
641fb90e4d1STodd Munson   tr->mu2_i = 0.50;
642fb90e4d1STodd Munson 
643fb90e4d1STodd Munson   tr->gamma1_i = 0.0625;
644fb90e4d1STodd Munson   tr->gamma2_i = 0.50;
645fb90e4d1STodd Munson   tr->gamma3_i = 2.00;
646fb90e4d1STodd Munson   tr->gamma4_i = 5.00;
647fb90e4d1STodd Munson 
648fb90e4d1STodd Munson   tr->theta_i = 0.25;
649fb90e4d1STodd Munson 
650a7e14dcfSSatish Balay   tr->min_radius = 1.0e-10;
651a7e14dcfSSatish Balay   tr->max_radius = 1.0e10;
652a7e14dcfSSatish Balay   tr->epsilon    = 1.0e-6;
653a7e14dcfSSatish Balay 
654a7e14dcfSSatish Balay   tr->init_type       = NTR_INIT_INTERPOLATION;
655a7e14dcfSSatish Balay   tr->update_type     = NTR_UPDATE_REDUCTION;
656a7e14dcfSSatish Balay 
657a7e14dcfSSatish Balay   /* Set linear solver to default for trust region */
658a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
65963b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1);CHKERRQ(ierr);
6605d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
661cbf034f8SAlp Dener   ierr = KSPAppendOptionsPrefix(tao->ksp,"tao_ntr_");CHKERRQ(ierr);
66205de396fSBarry Smith   ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr);
663a7e14dcfSSatish Balay   PetscFunctionReturn(0);
664a7e14dcfSSatish Balay }
665