xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision cbf034f8d8e89e62e6f33df37e349a71dad3c0f3)
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 
39ba7fe8fbSTodd Munson    Note:  TaoSolve_NTR MUST use the iterative solver KSPCGNASH, KSPCGSTCG,
40ba7fe8fbSTodd Munson           or KSPCGGLTR.  Thus, we set KSPCGNASH, KSPCGSTCG, or KSPCGGLTR 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) {
65a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"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);
69fb90e4d1STodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&is_nash);CHKERRQ(ierr);
70fb90e4d1STodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&is_stcg);CHKERRQ(ierr);
71fb90e4d1STodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&is_gltr);CHKERRQ(ierr);
72fb90e4d1STodd Munson   if (!is_nash && !is_stcg && !is_gltr) {
73fb90e4d1STodd Munson     SETERRQ(PETSC_COMM_SELF,1,"TAO_NTR requires nash, stcg, or gltr for the KSP");
74fb90e4d1STodd Munson   }
75a7e14dcfSSatish Balay 
76fb90e4d1STodd Munson   /* Initialize the radius and modify if it is too large or small */
77fb90e4d1STodd Munson   tao->trust = tao->trust0;
78a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tr->min_radius);
79a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tr->max_radius);
80a7e14dcfSSatish Balay 
810c51296cSAlp Dener /* Allocate the vectors needed for the BFGS approximation */
820c51296cSAlp Dener ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
830c51296cSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr);
840c51296cSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr);
850c51296cSAlp Dener if (is_bfgs) {
860c51296cSAlp Dener   tr->bfgs_pre = pc;
870c51296cSAlp Dener   ierr = PCLMVMGetMatLMVM(tr->bfgs_pre, &tr->M);CHKERRQ(ierr);
88a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
89a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
900c51296cSAlp Dener   ierr = MatSetSizes(tr->M, n, n, N, N);CHKERRQ(ierr);
91cd929ea3SAlp Dener   ierr = MatLMVMAllocate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
920ad3a497SAlp Dener   ierr = MatIsSymmetricKnown(tr->M, &sym_set, &is_symmetric);CHKERRQ(ierr);
930ad3a497SAlp Dener   if (!sym_set || !is_symmetric) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
940c51296cSAlp Dener } else if (is_jacobi) {
950c51296cSAlp Dener   ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
96a7e14dcfSSatish Balay }
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   /* Check convergence criteria */
99a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
100a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
10153506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Inf or NaN");
102a7e14dcfSSatish Balay   needH = 1;
103a7e14dcfSSatish Balay 
1043ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1053ecd9318SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1063ecd9318SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,1.0);CHKERRQ(ierr);
1073ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1083ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay   /*  Initialize trust-region radius */
111a7e14dcfSSatish Balay   switch(tr->init_type) {
112a7e14dcfSSatish Balay   case NTR_INIT_CONSTANT:
113a7e14dcfSSatish Balay     /*  Use the initial radius specified */
114a7e14dcfSSatish Balay     break;
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay   case NTR_INIT_INTERPOLATION:
117a7e14dcfSSatish Balay     /*  Use the initial radius specified */
118a7e14dcfSSatish Balay     max_radius = 0.0;
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
121a7e14dcfSSatish Balay       fmin = f;
122a7e14dcfSSatish Balay       sigma = 0.0;
123a7e14dcfSSatish Balay 
124a7e14dcfSSatish Balay       if (needH) {
125ffad9901SBarry Smith         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
126a7e14dcfSSatish Balay         needH = 0;
127a7e14dcfSSatish Balay       }
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
132a7e14dcfSSatish Balay         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
133a7e14dcfSSatish Balay         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
136a7e14dcfSSatish Balay           tau = tr->gamma1_i;
137a7e14dcfSSatish Balay         }
138a7e14dcfSSatish Balay         else {
139a7e14dcfSSatish Balay           if (ftrial < fmin) {
140a7e14dcfSSatish Balay             fmin = ftrial;
141a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
142a7e14dcfSSatish Balay           }
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
145a7e14dcfSSatish Balay           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
148a7e14dcfSSatish Balay           actred = f - ftrial;
149a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
150a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tr->epsilon)) {
151a7e14dcfSSatish Balay             kappa = 1.0;
152a7e14dcfSSatish Balay           }
153a7e14dcfSSatish Balay           else {
154a7e14dcfSSatish Balay             kappa = actred / prered;
155a7e14dcfSSatish Balay           }
156a7e14dcfSSatish Balay 
157a7e14dcfSSatish Balay           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
158a7e14dcfSSatish Balay           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
159a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
160a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
161a7e14dcfSSatish Balay 
162a7e14dcfSSatish Balay           if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) {
163a7e14dcfSSatish Balay             /*  Great agreement */
164a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
165a7e14dcfSSatish Balay 
166a7e14dcfSSatish Balay             if (tau_max < 1.0) {
167a7e14dcfSSatish Balay               tau = tr->gamma3_i;
168a7e14dcfSSatish Balay             }
169a7e14dcfSSatish Balay             else if (tau_max > tr->gamma4_i) {
170a7e14dcfSSatish Balay               tau = tr->gamma4_i;
171a7e14dcfSSatish Balay             }
172a7e14dcfSSatish Balay             else {
173a7e14dcfSSatish Balay               tau = tau_max;
174a7e14dcfSSatish Balay             }
175a7e14dcfSSatish Balay           }
176a7e14dcfSSatish Balay           else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) {
177a7e14dcfSSatish Balay             /*  Good agreement */
178a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
179a7e14dcfSSatish Balay 
180a7e14dcfSSatish Balay             if (tau_max < tr->gamma2_i) {
181a7e14dcfSSatish Balay               tau = tr->gamma2_i;
182a7e14dcfSSatish Balay             }
183a7e14dcfSSatish Balay             else if (tau_max > tr->gamma3_i) {
184a7e14dcfSSatish Balay               tau = tr->gamma3_i;
185a7e14dcfSSatish Balay             }
186a7e14dcfSSatish Balay             else {
187a7e14dcfSSatish Balay               tau = tau_max;
188a7e14dcfSSatish Balay             }
189a7e14dcfSSatish Balay           }
190a7e14dcfSSatish Balay           else {
191a7e14dcfSSatish Balay             /*  Not good agreement */
192a7e14dcfSSatish Balay             if (tau_min > 1.0) {
193a7e14dcfSSatish Balay               tau = tr->gamma2_i;
194a7e14dcfSSatish Balay             }
195a7e14dcfSSatish Balay             else if (tau_max < tr->gamma1_i) {
196a7e14dcfSSatish Balay               tau = tr->gamma1_i;
197a7e14dcfSSatish Balay             }
198a7e14dcfSSatish Balay             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
199a7e14dcfSSatish Balay               tau = tr->gamma1_i;
200a7e14dcfSSatish Balay             }
201a7e14dcfSSatish Balay             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
202a7e14dcfSSatish Balay                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
203a7e14dcfSSatish Balay               tau = tau_1;
204a7e14dcfSSatish Balay             }
205a7e14dcfSSatish Balay             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
206a7e14dcfSSatish Balay                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
207a7e14dcfSSatish Balay               tau = tau_2;
208a7e14dcfSSatish Balay             }
209a7e14dcfSSatish Balay             else {
210a7e14dcfSSatish Balay               tau = tau_max;
211a7e14dcfSSatish Balay             }
212a7e14dcfSSatish Balay           }
213a7e14dcfSSatish Balay         }
214a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
215a7e14dcfSSatish Balay       }
216a7e14dcfSSatish Balay 
217a7e14dcfSSatish Balay       if (fmin < f) {
218a7e14dcfSSatish Balay         f = fmin;
219a7e14dcfSSatish Balay         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
220a7e14dcfSSatish Balay         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
221a7e14dcfSSatish Balay 
222a9603a14SPatrick Farrell         ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
223a7e14dcfSSatish Balay 
22453506e15SBarry Smith         if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
225a7e14dcfSSatish Balay         needH = 1;
226a7e14dcfSSatish Balay 
2273ecd9318SAlp Dener         ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
2283ecd9318SAlp Dener         ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,1.0);CHKERRQ(ierr);
2293ecd9318SAlp Dener         ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2303ecd9318SAlp Dener         if (tao->reason != TAO_CONTINUE_ITERATING) {
231a7e14dcfSSatish Balay           PetscFunctionReturn(0);
232a7e14dcfSSatish Balay         }
233a7e14dcfSSatish Balay       }
234a7e14dcfSSatish Balay     }
235a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay     /*  Modify the radius if it is too large or small */
238a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tr->min_radius);
239a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
240a7e14dcfSSatish Balay     break;
241a7e14dcfSSatish Balay 
242a7e14dcfSSatish Balay   default:
243a7e14dcfSSatish Balay     /*  Norm of the first direction will initialize radius */
244a7e14dcfSSatish Balay     tao->trust = 0.0;
245a7e14dcfSSatish Balay     break;
246a7e14dcfSSatish Balay   }
247a7e14dcfSSatish Balay 
248a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2493ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
2508931d482SJason Sarich     ++tao->niter;
251ae93cb3cSJason Sarich     tao->ksp_its=0;
252a7e14dcfSSatish Balay     /* Compute the Hessian */
253a7e14dcfSSatish Balay     if (needH) {
254ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
255a7e14dcfSSatish Balay       needH = 0;
256a7e14dcfSSatish Balay     }
257a7e14dcfSSatish Balay 
2580c51296cSAlp Dener     if (tr->bfgs_pre) {
259a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
260a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
261a7e14dcfSSatish Balay       ++bfgsUpdates;
262a7e14dcfSSatish Balay     }
263a7e14dcfSSatish Balay 
2643ecd9318SAlp Dener     while (tao->reason == TAO_CONTINUE_ITERATING) {
26523ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
266a7e14dcfSSatish Balay 
267a7e14dcfSSatish Balay       /* Solve the trust region subproblem */
268ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
269a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
270a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
271a7e14dcfSSatish Balay       tao->ksp_its+=its;
272ae93cb3cSJason Sarich       tao->ksp_tot_its+=its;
273ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
274a7e14dcfSSatish Balay 
275a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
276a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
277a7e14dcfSSatish Balay         if (norm_d > 0.0) {
278a7e14dcfSSatish Balay           tao->trust = norm_d;
279a7e14dcfSSatish Balay 
280a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
281a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
282a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
283a7e14dcfSSatish Balay         }
284a7e14dcfSSatish Balay         else {
285a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
286a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
287a7e14dcfSSatish Balay           tao->trust = tao->trust0;
288a7e14dcfSSatish Balay 
289a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
290a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
291a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
292a7e14dcfSSatish Balay 
293ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
294a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
295a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
296a7e14dcfSSatish Balay           tao->ksp_its+=its;
2972d9aa51bSJason Sarich           tao->ksp_tot_its+=its;
298ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
299a7e14dcfSSatish Balay 
30053506e15SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
301a7e14dcfSSatish Balay         }
302a7e14dcfSSatish Balay       }
303a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
304a7e14dcfSSatish Balay       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
3050c51296cSAlp Dener       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (tr->bfgs_pre)) {
306a7e14dcfSSatish Balay         /* Preconditioner is numerically indefinite; reset the
307a7e14dcfSSatish Balay            approximate if using BFGS preconditioning. */
308cd929ea3SAlp Dener         ierr = MatLMVMReset(tr->M, PETSC_FALSE);CHKERRQ(ierr);
309a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
310a7e14dcfSSatish Balay         bfgsUpdates = 1;
311a7e14dcfSSatish Balay       }
312a7e14dcfSSatish Balay 
313a7e14dcfSSatish Balay       if (NTR_UPDATE_REDUCTION == tr->update_type) {
314a7e14dcfSSatish Balay         /* Get predicted reduction */
315ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
316a7e14dcfSSatish Balay         if (prered >= 0.0) {
317a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
318a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
319a7e14dcfSSatish Balay              be rejected! */
320a7e14dcfSSatish Balay           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
321a7e14dcfSSatish Balay         }
322a7e14dcfSSatish Balay         else {
323a7e14dcfSSatish Balay           /* Compute trial step and function value */
324a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
325a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
326a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
327a7e14dcfSSatish Balay 
328a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
329a7e14dcfSSatish Balay             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
330a7e14dcfSSatish Balay           } else {
331a7e14dcfSSatish Balay             /* Compute and actual reduction */
332a7e14dcfSSatish Balay             actred = f - ftrial;
333a7e14dcfSSatish Balay             prered = -prered;
334a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
335a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
336a7e14dcfSSatish Balay               kappa = 1.0;
337a7e14dcfSSatish Balay             }
338a7e14dcfSSatish Balay             else {
339a7e14dcfSSatish Balay               kappa = actred / prered;
340a7e14dcfSSatish Balay             }
341a7e14dcfSSatish Balay 
342a7e14dcfSSatish Balay             /* Accept or reject the step and update radius */
343a7e14dcfSSatish Balay             if (kappa < tr->eta1) {
344a7e14dcfSSatish Balay               /* Reject the step */
345a7e14dcfSSatish Balay               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
346a7e14dcfSSatish Balay             }
347a7e14dcfSSatish Balay             else {
348a7e14dcfSSatish Balay               /* Accept the step */
349a7e14dcfSSatish Balay               if (kappa < tr->eta2) {
350a7e14dcfSSatish Balay                 /* Marginal bad step */
351a7e14dcfSSatish Balay                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
352a7e14dcfSSatish Balay               }
353a7e14dcfSSatish Balay               else if (kappa < tr->eta3) {
354a7e14dcfSSatish Balay                 /* Reasonable step */
355a7e14dcfSSatish Balay                 tao->trust = tr->alpha3 * tao->trust;
356a7e14dcfSSatish Balay               }
357a7e14dcfSSatish Balay               else if (kappa < tr->eta4) {
358a7e14dcfSSatish Balay                 /* Good step */
359a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
360a7e14dcfSSatish Balay               }
361a7e14dcfSSatish Balay               else {
362a7e14dcfSSatish Balay                 /* Very good step */
363a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
364a7e14dcfSSatish Balay               }
365a7e14dcfSSatish Balay               break;
366a7e14dcfSSatish Balay             }
367a7e14dcfSSatish Balay           }
368a7e14dcfSSatish Balay         }
369a7e14dcfSSatish Balay       }
370a7e14dcfSSatish Balay       else {
371a7e14dcfSSatish Balay         /* Get predicted reduction */
372ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
373a7e14dcfSSatish Balay         if (prered >= 0.0) {
374a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
375a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
376a7e14dcfSSatish Balay              be rejected! */
377a7e14dcfSSatish Balay           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
378a7e14dcfSSatish Balay         }
379a7e14dcfSSatish Balay         else {
380a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
381a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
382a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
383a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
384a7e14dcfSSatish Balay             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
385a7e14dcfSSatish Balay           }
386a7e14dcfSSatish Balay           else {
387a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
388a7e14dcfSSatish Balay             actred = f - ftrial;
389a7e14dcfSSatish Balay             prered = -prered;
390a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
391a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
392a7e14dcfSSatish Balay               kappa = 1.0;
393a7e14dcfSSatish Balay             }
394a7e14dcfSSatish Balay             else {
395a7e14dcfSSatish Balay               kappa = actred / prered;
396a7e14dcfSSatish Balay             }
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
399a7e14dcfSSatish Balay             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
400a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
401a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
402a7e14dcfSSatish Balay 
403a7e14dcfSSatish Balay             if (kappa >= 1.0 - tr->mu1) {
404a7e14dcfSSatish Balay               /* Great agreement; accept step and update radius */
405a7e14dcfSSatish Balay               if (tau_max < 1.0) {
406a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
407a7e14dcfSSatish Balay               }
408a7e14dcfSSatish Balay               else if (tau_max > tr->gamma4) {
409a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
410a7e14dcfSSatish Balay               }
411a7e14dcfSSatish Balay               else {
412a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
413a7e14dcfSSatish Balay               }
414a7e14dcfSSatish Balay               break;
415a7e14dcfSSatish Balay             }
416a7e14dcfSSatish Balay             else if (kappa >= 1.0 - tr->mu2) {
417a7e14dcfSSatish Balay               /* Good agreement */
418a7e14dcfSSatish Balay 
419a7e14dcfSSatish Balay               if (tau_max < tr->gamma2) {
420a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
421a7e14dcfSSatish Balay               }
422a7e14dcfSSatish Balay               else if (tau_max > tr->gamma3) {
423a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
424a7e14dcfSSatish Balay               }
425a7e14dcfSSatish Balay               else if (tau_max < 1.0) {
426a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
427a7e14dcfSSatish Balay               }
428a7e14dcfSSatish Balay               else {
429a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
430a7e14dcfSSatish Balay               }
431a7e14dcfSSatish Balay               break;
432a7e14dcfSSatish Balay             }
433a7e14dcfSSatish Balay             else {
434a7e14dcfSSatish Balay               /* Not good agreement */
435a7e14dcfSSatish Balay               if (tau_min > 1.0) {
436a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
437a7e14dcfSSatish Balay               }
438a7e14dcfSSatish Balay               else if (tau_max < tr->gamma1) {
439a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
440a7e14dcfSSatish Balay               }
441a7e14dcfSSatish Balay               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
442a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
443a7e14dcfSSatish Balay               }
444a7e14dcfSSatish Balay               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
445a7e14dcfSSatish Balay                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
446a7e14dcfSSatish Balay                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
447a7e14dcfSSatish Balay               }
448a7e14dcfSSatish Balay               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
449a7e14dcfSSatish Balay                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
450a7e14dcfSSatish Balay                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
451a7e14dcfSSatish Balay               }
452a7e14dcfSSatish Balay               else {
453a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
454a7e14dcfSSatish Balay               }
455a7e14dcfSSatish Balay             }
456a7e14dcfSSatish Balay           }
457a7e14dcfSSatish Balay         }
458a7e14dcfSSatish Balay       }
459a7e14dcfSSatish Balay 
460a7e14dcfSSatish Balay       /* The step computed was not good and the radius was decreased.
461a7e14dcfSSatish Balay          Monitor the radius to terminate. */
4623ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
4633ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,tao->trust);CHKERRQ(ierr);
4643ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
465a7e14dcfSSatish Balay     }
466a7e14dcfSSatish Balay 
467a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
468a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
469a7e14dcfSSatish Balay 
4703ecd9318SAlp Dener     if (tao->reason == TAO_CONTINUE_ITERATING) {
471a7e14dcfSSatish Balay       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
472a7e14dcfSSatish Balay       f = ftrial;
473302440fdSBarry Smith       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr);
474a9603a14SPatrick Farrell       ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
47553506e15SBarry Smith       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
476a7e14dcfSSatish Balay       needH = 1;
4773ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
4783ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,tao->trust);CHKERRQ(ierr);
4793ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
480a7e14dcfSSatish Balay     }
481a7e14dcfSSatish Balay   }
482a7e14dcfSSatish Balay   PetscFunctionReturn(0);
483a7e14dcfSSatish Balay }
484a7e14dcfSSatish Balay 
485a7e14dcfSSatish Balay /*------------------------------------------------------------*/
486441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao)
487a7e14dcfSSatish Balay {
488a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
489a7e14dcfSSatish Balay   PetscErrorCode ierr;
490a7e14dcfSSatish Balay 
491a7e14dcfSSatish Balay   PetscFunctionBegin;
492a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
493a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
494a7e14dcfSSatish Balay   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
495a7e14dcfSSatish Balay 
4960c51296cSAlp Dener   tr->bfgs_pre = 0;
497a7e14dcfSSatish Balay   tr->M = 0;
498a7e14dcfSSatish Balay   PetscFunctionReturn(0);
499a7e14dcfSSatish Balay }
500a7e14dcfSSatish Balay 
501a7e14dcfSSatish Balay /*------------------------------------------------------------*/
502441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTR(Tao tao)
503a7e14dcfSSatish Balay {
504a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
505a7e14dcfSSatish Balay   PetscErrorCode ierr;
506a7e14dcfSSatish Balay 
507a7e14dcfSSatish Balay   PetscFunctionBegin;
508a7e14dcfSSatish Balay   if (tao->setupcalled) {
509a7e14dcfSSatish Balay     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
510a7e14dcfSSatish Balay   }
511a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
512a7e14dcfSSatish Balay   PetscFunctionReturn(0);
513a7e14dcfSSatish Balay }
514a7e14dcfSSatish Balay 
515a7e14dcfSSatish Balay /*------------------------------------------------------------*/
5164416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(PetscOptionItems *PetscOptionsObject,Tao tao)
517a7e14dcfSSatish Balay {
518a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
519a7e14dcfSSatish Balay   PetscErrorCode ierr;
520a7e14dcfSSatish Balay 
521a7e14dcfSSatish Balay   PetscFunctionBegin;
5221a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
52394ae4db5SBarry 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);
52494ae4db5SBarry 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);
52594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr);
52694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr);
52794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr);
52894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr);
52994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr);
53094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr);
53194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr);
53294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr);
53394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr);
53494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr);
53594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr);
53694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr);
53794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr);
53894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr);
53994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr);
54094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr);
54194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr);
54294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr);
54394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr);
54494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr);
54594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr);
54694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr);
54794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr);
54894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr);
54994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr);
55094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr);
551a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
552a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
553a7e14dcfSSatish Balay   PetscFunctionReturn(0);
554a7e14dcfSSatish Balay }
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay /*------------------------------------------------------------*/
5571522df2eSJason Sarich /*MC
5581522df2eSJason Sarich   TAONTR - Newton's method with trust region for unconstrained minimization.
5591522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system.
5601522df2eSJason Sarich   NTR expects a KSP solver with a trust region radius.
5611522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
5621522df2eSJason Sarich 
5631522df2eSJason Sarich   Options Database Keys:
564fb90e4d1STodd Munson + -tao_ntr_pc_type - "none","ahess","bfgs","petsc"
5651522df2eSJason Sarich . -tao_ntr_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs"
5661522df2eSJason Sarich . -tao_ntr_init_type - "constant","direction","interpolation"
5671522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation"
5681522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius
5691522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius
5701522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
5711522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor
5721522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor
5731522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor
5741522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor
5751522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor
5761522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor
5771522df2eSJason Sarich . -tao_ntr_theta_i - thetha1 interpolation init factor
5781522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor
5791522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor
5801522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor
5811522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor
5821522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor
5831522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor
5841522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor
5851522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
5861522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
5871522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update
5881522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update
5891522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update
5901522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update
5911522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update
5921522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update
5931522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update
5941522df2eSJason Sarich 
5951eb8069cSJason Sarich   Level: beginner
5961522df2eSJason Sarich M*/
5971522df2eSJason Sarich 
598728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
599a7e14dcfSSatish Balay {
600a7e14dcfSSatish Balay   TAO_NTR *tr;
601a7e14dcfSSatish Balay   PetscErrorCode ierr;
602a7e14dcfSSatish Balay 
603a7e14dcfSSatish Balay   PetscFunctionBegin;
604a7e14dcfSSatish Balay 
6053c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
606a7e14dcfSSatish Balay 
607a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTR;
608a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTR;
609a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
610a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTR;
611a7e14dcfSSatish Balay 
6126552cf8aSJason Sarich   /* Override default settings (unless already changed) */
6136552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
6146552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
615a7e14dcfSSatish Balay   tao->data = (void*)tr;
616a7e14dcfSSatish Balay 
617a7e14dcfSSatish Balay   /*  Standard trust region update parameters */
618a7e14dcfSSatish Balay   tr->eta1 = 1.0e-4;
619a7e14dcfSSatish Balay   tr->eta2 = 0.25;
620a7e14dcfSSatish Balay   tr->eta3 = 0.50;
621a7e14dcfSSatish Balay   tr->eta4 = 0.90;
622a7e14dcfSSatish Balay 
623a7e14dcfSSatish Balay   tr->alpha1 = 0.25;
624a7e14dcfSSatish Balay   tr->alpha2 = 0.50;
625a7e14dcfSSatish Balay   tr->alpha3 = 1.00;
626a7e14dcfSSatish Balay   tr->alpha4 = 2.00;
627a7e14dcfSSatish Balay   tr->alpha5 = 4.00;
628a7e14dcfSSatish Balay 
629a7e14dcfSSatish Balay   /*  Interpolation trust region update parameters */
630a7e14dcfSSatish Balay   tr->mu1 = 0.10;
631a7e14dcfSSatish Balay   tr->mu2 = 0.50;
632a7e14dcfSSatish Balay 
633a7e14dcfSSatish Balay   tr->gamma1 = 0.25;
634a7e14dcfSSatish Balay   tr->gamma2 = 0.50;
635a7e14dcfSSatish Balay   tr->gamma3 = 2.00;
636a7e14dcfSSatish Balay   tr->gamma4 = 4.00;
637a7e14dcfSSatish Balay 
638a7e14dcfSSatish Balay   tr->theta = 0.05;
639a7e14dcfSSatish Balay 
640fb90e4d1STodd Munson   /*  Interpolation parameters for initialization */
641fb90e4d1STodd Munson   tr->mu1_i = 0.35;
642fb90e4d1STodd Munson   tr->mu2_i = 0.50;
643fb90e4d1STodd Munson 
644fb90e4d1STodd Munson   tr->gamma1_i = 0.0625;
645fb90e4d1STodd Munson   tr->gamma2_i = 0.50;
646fb90e4d1STodd Munson   tr->gamma3_i = 2.00;
647fb90e4d1STodd Munson   tr->gamma4_i = 5.00;
648fb90e4d1STodd Munson 
649fb90e4d1STodd Munson   tr->theta_i = 0.25;
650fb90e4d1STodd Munson 
651a7e14dcfSSatish Balay   tr->min_radius = 1.0e-10;
652a7e14dcfSSatish Balay   tr->max_radius = 1.0e10;
653a7e14dcfSSatish Balay   tr->epsilon    = 1.0e-6;
654a7e14dcfSSatish Balay 
655a7e14dcfSSatish Balay   tr->init_type       = NTR_INIT_INTERPOLATION;
656a7e14dcfSSatish Balay   tr->update_type     = NTR_UPDATE_REDUCTION;
657a7e14dcfSSatish Balay 
658a7e14dcfSSatish Balay   /* Set linear solver to default for trust region */
659a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
66063b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
661*cbf034f8SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
662*cbf034f8SAlp Dener   ierr = KSPAppendOptionsPrefix(tao->ksp, "tao_ntr_");CHKERRQ(ierr);
663fb90e4d1STodd Munson   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
664a7e14dcfSSatish Balay   PetscFunctionReturn(0);
665a7e14dcfSSatish Balay }
666