xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision cd929ea3f739fd9f7b6394f772cb40b9d4e6d97c)
1fb90e4d1STodd Munson #include <../src/tao/unconstrained/impls/ntr/ntrimpl.h>
2a7e14dcfSSatish Balay 
3aaa7dc30SBarry Smith #include <petscksp.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay #define NTR_PC_NONE     0
6a7e14dcfSSatish Balay #define NTR_PC_AHESS    1
7a7e14dcfSSatish Balay #define NTR_PC_BFGS     2
8a7e14dcfSSatish Balay #define NTR_PC_PETSC    3
9a7e14dcfSSatish Balay #define NTR_PC_TYPES    4
10a7e14dcfSSatish Balay 
11a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS   0
12a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS    1
13a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES   2
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay #define NTR_INIT_CONSTANT         0
16a7e14dcfSSatish Balay #define NTR_INIT_DIRECTION        1
17a7e14dcfSSatish Balay #define NTR_INIT_INTERPOLATION    2
18a7e14dcfSSatish Balay #define NTR_INIT_TYPES            3
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay #define NTR_UPDATE_REDUCTION      0
21a7e14dcfSSatish Balay #define NTR_UPDATE_INTERPOLATION  1
22a7e14dcfSSatish Balay #define NTR_UPDATE_TYPES          2
23a7e14dcfSSatish Balay 
2453506e15SBarry Smith static const char *NTR_PC[64] = {"none","ahess","bfgs","petsc"};
25a7e14dcfSSatish Balay 
2653506e15SBarry Smith static const char *BFGS_SCALE[64] = {"ahess","bfgs"};
27a7e14dcfSSatish Balay 
2853506e15SBarry Smith static const char *NTR_INIT[64] = {"constant","direction","interpolation"};
29a7e14dcfSSatish Balay 
3053506e15SBarry Smith static const char *NTR_UPDATE[64] = {"reduction","interpolation"};
31a7e14dcfSSatish Balay 
32*cd929ea3SAlp Dener PetscErrorCode TaoNTRPreconBFGS(PC BFGSpc, Vec X, Vec Y)
33fb90e4d1STodd Munson {
34fb90e4d1STodd Munson   PetscErrorCode ierr;
35*cd929ea3SAlp Dener   Mat *M;
36*cd929ea3SAlp Dener 
37fb90e4d1STodd Munson   PetscFunctionBegin;
38*cd929ea3SAlp Dener   ierr = PCShellGetContext(BFGSpc, (void**)&M);CHKERRQ(ierr);
39*cd929ea3SAlp Dener   ierr = MatSolve(*M, X, Y);CHKERRQ(ierr);
40fb90e4d1STodd Munson   PetscFunctionReturn(0);
41fb90e4d1STodd Munson }
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 
64ba7fe8fbSTodd Munson    Note:  TaoSolve_NTR MUST use the iterative solver KSPCGNASH, KSPCGSTCG,
65ba7fe8fbSTodd 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 */
68441846f8SBarry Smith static PetscErrorCode TaoSolve_NTR(Tao tao)
69a7e14dcfSSatish Balay {
70a7e14dcfSSatish Balay   TAO_NTR            *tr = (TAO_NTR *)tao->data;
71fb90e4d1STodd Munson   KSPType            ksp_type;
72fb90e4d1STodd Munson   PetscBool          is_nash,is_stcg,is_gltr;
73a7e14dcfSSatish Balay   KSPConvergedReason ksp_reason;
74fb90e4d1STodd Munson   PC                 pc;
75a7e14dcfSSatish Balay   PetscReal          fmin, ftrial, prered, actred, kappa, sigma, beta;
76a7e14dcfSSatish Balay   PetscReal          tau, tau_1, tau_2, tau_max, tau_min, max_radius;
77a7e14dcfSSatish Balay   PetscReal          f, gnorm;
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay   PetscReal          delta;
80a7e14dcfSSatish Balay   PetscReal          norm_d;
81a7e14dcfSSatish Balay   PetscErrorCode     ierr;
82a7e14dcfSSatish Balay   PetscInt           bfgsUpdates = 0;
83a7e14dcfSSatish Balay   PetscInt           needH;
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay   PetscInt           i_max = 5;
86a7e14dcfSSatish Balay   PetscInt           j_max = 1;
87a7e14dcfSSatish Balay   PetscInt           i, j, N, n, its;
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay   PetscFunctionBegin;
90a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
91a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr);
92a7e14dcfSSatish Balay   }
93a7e14dcfSSatish Balay 
94fb90e4d1STodd Munson   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
95fb90e4d1STodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&is_nash);CHKERRQ(ierr);
96fb90e4d1STodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&is_stcg);CHKERRQ(ierr);
97fb90e4d1STodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&is_gltr);CHKERRQ(ierr);
98fb90e4d1STodd Munson   if (!is_nash && !is_stcg && !is_gltr) {
99fb90e4d1STodd Munson     SETERRQ(PETSC_COMM_SELF,1,"TAO_NTR requires nash, stcg, or gltr for the KSP");
100fb90e4d1STodd Munson   }
101a7e14dcfSSatish Balay 
102fb90e4d1STodd Munson   /* Initialize the radius and modify if it is too large or small */
103fb90e4d1STodd Munson   tao->trust = tao->trust0;
104a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tr->min_radius);
105a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tr->max_radius);
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type && !tr->M) {
108a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
109a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
110*cd929ea3SAlp Dener     ierr = MatCreateLBFGS(((PetscObject)tao)->comm,n,N,&tr->M);CHKERRQ(ierr);
111*cd929ea3SAlp Dener     ierr = MatLMVMAllocate(tr->M,tao->solution,tao->gradient);CHKERRQ(ierr);
112a7e14dcfSSatish Balay   }
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay   /* Check convergence criteria */
115a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
116a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
11753506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Inf or NaN");
118a7e14dcfSSatish Balay   needH = 1;
119a7e14dcfSSatish Balay 
1203ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1213ecd9318SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1223ecd9318SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,1.0);CHKERRQ(ierr);
1233ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1243ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay   /* Create vectors for the limited memory preconditioner */
127fb90e4d1STodd Munson   if ((NTR_PC_BFGS == tr->pc_type) && (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) {
128a7e14dcfSSatish Balay     if (!tr->Diag) {
129a7e14dcfSSatish Balay         ierr = VecDuplicate(tao->solution, &tr->Diag);CHKERRQ(ierr);
130a7e14dcfSSatish Balay     }
131a7e14dcfSSatish Balay   }
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay   /*  Modify the preconditioner to use the bfgs approximation */
134a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
135a7e14dcfSSatish Balay   switch(tr->pc_type) {
136a7e14dcfSSatish Balay   case NTR_PC_NONE:
137a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1381a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
139a7e14dcfSSatish Balay     break;
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay   case NTR_PC_AHESS:
142a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1431a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
144baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
145a7e14dcfSSatish Balay     break;
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay   case NTR_PC_BFGS:
148a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
1491a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
150a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
151a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, tr->M);CHKERRQ(ierr);
152*cd929ea3SAlp Dener     ierr = PCShellSetApply(pc, TaoNTRPreconBFGS);CHKERRQ(ierr);
153a7e14dcfSSatish Balay     break;
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay   default:
156a7e14dcfSSatish Balay     /*  Use the pc method set by pc_type */
157a7e14dcfSSatish Balay     break;
158a7e14dcfSSatish Balay   }
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay   /*  Initialize trust-region radius */
161a7e14dcfSSatish Balay   switch(tr->init_type) {
162a7e14dcfSSatish Balay   case NTR_INIT_CONSTANT:
163a7e14dcfSSatish Balay     /*  Use the initial radius specified */
164a7e14dcfSSatish Balay     break;
165a7e14dcfSSatish Balay 
166a7e14dcfSSatish Balay   case NTR_INIT_INTERPOLATION:
167a7e14dcfSSatish Balay     /*  Use the initial radius specified */
168a7e14dcfSSatish Balay     max_radius = 0.0;
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
171a7e14dcfSSatish Balay       fmin = f;
172a7e14dcfSSatish Balay       sigma = 0.0;
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay       if (needH) {
175ffad9901SBarry Smith         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
176a7e14dcfSSatish Balay         needH = 0;
177a7e14dcfSSatish Balay       }
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
182a7e14dcfSSatish Balay         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
183a7e14dcfSSatish Balay         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
186a7e14dcfSSatish Balay           tau = tr->gamma1_i;
187a7e14dcfSSatish Balay         }
188a7e14dcfSSatish Balay         else {
189a7e14dcfSSatish Balay           if (ftrial < fmin) {
190a7e14dcfSSatish Balay             fmin = ftrial;
191a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
192a7e14dcfSSatish Balay           }
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
195a7e14dcfSSatish Balay           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
196a7e14dcfSSatish Balay 
197a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
198a7e14dcfSSatish Balay           actred = f - ftrial;
199a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
200a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tr->epsilon)) {
201a7e14dcfSSatish Balay             kappa = 1.0;
202a7e14dcfSSatish Balay           }
203a7e14dcfSSatish Balay           else {
204a7e14dcfSSatish Balay             kappa = actred / prered;
205a7e14dcfSSatish Balay           }
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
208a7e14dcfSSatish Balay           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
209a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
210a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay           if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) {
213a7e14dcfSSatish Balay             /*  Great agreement */
214a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay             if (tau_max < 1.0) {
217a7e14dcfSSatish Balay               tau = tr->gamma3_i;
218a7e14dcfSSatish Balay             }
219a7e14dcfSSatish Balay             else if (tau_max > tr->gamma4_i) {
220a7e14dcfSSatish Balay               tau = tr->gamma4_i;
221a7e14dcfSSatish Balay             }
222a7e14dcfSSatish Balay             else {
223a7e14dcfSSatish Balay               tau = tau_max;
224a7e14dcfSSatish Balay             }
225a7e14dcfSSatish Balay           }
226a7e14dcfSSatish Balay           else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) {
227a7e14dcfSSatish Balay             /*  Good agreement */
228a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay             if (tau_max < tr->gamma2_i) {
231a7e14dcfSSatish Balay               tau = tr->gamma2_i;
232a7e14dcfSSatish Balay             }
233a7e14dcfSSatish Balay             else if (tau_max > tr->gamma3_i) {
234a7e14dcfSSatish Balay               tau = tr->gamma3_i;
235a7e14dcfSSatish Balay             }
236a7e14dcfSSatish Balay             else {
237a7e14dcfSSatish Balay               tau = tau_max;
238a7e14dcfSSatish Balay             }
239a7e14dcfSSatish Balay           }
240a7e14dcfSSatish Balay           else {
241a7e14dcfSSatish Balay             /*  Not good agreement */
242a7e14dcfSSatish Balay             if (tau_min > 1.0) {
243a7e14dcfSSatish Balay               tau = tr->gamma2_i;
244a7e14dcfSSatish Balay             }
245a7e14dcfSSatish Balay             else if (tau_max < tr->gamma1_i) {
246a7e14dcfSSatish Balay               tau = tr->gamma1_i;
247a7e14dcfSSatish Balay             }
248a7e14dcfSSatish Balay             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
249a7e14dcfSSatish Balay               tau = tr->gamma1_i;
250a7e14dcfSSatish Balay             }
251a7e14dcfSSatish Balay             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
252a7e14dcfSSatish Balay                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
253a7e14dcfSSatish Balay               tau = tau_1;
254a7e14dcfSSatish Balay             }
255a7e14dcfSSatish Balay             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
256a7e14dcfSSatish Balay                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
257a7e14dcfSSatish Balay               tau = tau_2;
258a7e14dcfSSatish Balay             }
259a7e14dcfSSatish Balay             else {
260a7e14dcfSSatish Balay               tau = tau_max;
261a7e14dcfSSatish Balay             }
262a7e14dcfSSatish Balay           }
263a7e14dcfSSatish Balay         }
264a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
265a7e14dcfSSatish Balay       }
266a7e14dcfSSatish Balay 
267a7e14dcfSSatish Balay       if (fmin < f) {
268a7e14dcfSSatish Balay         f = fmin;
269a7e14dcfSSatish Balay         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
270a7e14dcfSSatish Balay         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
271a7e14dcfSSatish Balay 
272a9603a14SPatrick Farrell         ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
273a7e14dcfSSatish Balay 
27453506e15SBarry Smith         if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
275a7e14dcfSSatish Balay         needH = 1;
276a7e14dcfSSatish Balay 
2773ecd9318SAlp Dener         ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
2783ecd9318SAlp Dener         ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,1.0);CHKERRQ(ierr);
2793ecd9318SAlp Dener         ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2803ecd9318SAlp Dener         if (tao->reason != TAO_CONTINUE_ITERATING) {
281a7e14dcfSSatish Balay           PetscFunctionReturn(0);
282a7e14dcfSSatish Balay         }
283a7e14dcfSSatish Balay       }
284a7e14dcfSSatish Balay     }
285a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
286a7e14dcfSSatish Balay 
287a7e14dcfSSatish Balay     /*  Modify the radius if it is too large or small */
288a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tr->min_radius);
289a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
290a7e14dcfSSatish Balay     break;
291a7e14dcfSSatish Balay 
292a7e14dcfSSatish Balay   default:
293a7e14dcfSSatish Balay     /*  Norm of the first direction will initialize radius */
294a7e14dcfSSatish Balay     tao->trust = 0.0;
295a7e14dcfSSatish Balay     break;
296a7e14dcfSSatish Balay   }
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
299a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
300a7e14dcfSSatish Balay      since the function value may have decreased */
301a7e14dcfSSatish Balay   if (NTR_PC_BFGS == tr->pc_type) {
302*cd929ea3SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
303*cd929ea3SAlp Dener     ierr = MatLMVMSetJ0Scale(tr->M, delta);CHKERRQ(ierr);
304a7e14dcfSSatish Balay   }
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
3073ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
3088931d482SJason Sarich     ++tao->niter;
309ae93cb3cSJason Sarich     tao->ksp_its=0;
310a7e14dcfSSatish Balay     /* Compute the Hessian */
311a7e14dcfSSatish Balay     if (needH) {
312ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
313a7e14dcfSSatish Balay       needH = 0;
314a7e14dcfSSatish Balay     }
315a7e14dcfSSatish Balay 
316a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type) {
317a7e14dcfSSatish Balay       if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) {
318a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner */
319a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, tr->Diag);CHKERRQ(ierr);
320a7e14dcfSSatish Balay         ierr = VecAbs(tr->Diag);CHKERRQ(ierr);
321a7e14dcfSSatish Balay         ierr = VecReciprocal(tr->Diag);CHKERRQ(ierr);
322*cd929ea3SAlp Dener         ierr = MatLMVMSetJ0Diag(tr->M,tr->Diag);CHKERRQ(ierr);
323a7e14dcfSSatish Balay       }
324a7e14dcfSSatish Balay 
325a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
326a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
327a7e14dcfSSatish Balay       ++bfgsUpdates;
328a7e14dcfSSatish Balay     }
329a7e14dcfSSatish Balay 
3303ecd9318SAlp Dener     while (tao->reason == TAO_CONTINUE_ITERATING) {
33123ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
332a7e14dcfSSatish Balay 
333a7e14dcfSSatish Balay       /* Solve the trust region subproblem */
334ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
335a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
336a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
337a7e14dcfSSatish Balay       tao->ksp_its+=its;
338ae93cb3cSJason Sarich       tao->ksp_tot_its+=its;
339ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
340a7e14dcfSSatish Balay 
341a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
342a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
343a7e14dcfSSatish Balay         if (norm_d > 0.0) {
344a7e14dcfSSatish Balay           tao->trust = norm_d;
345a7e14dcfSSatish Balay 
346a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
347a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
348a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
349a7e14dcfSSatish Balay         }
350a7e14dcfSSatish Balay         else {
351a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
352a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
353a7e14dcfSSatish Balay           tao->trust = tao->trust0;
354a7e14dcfSSatish Balay 
355a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
356a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, tr->min_radius);
357a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, tr->max_radius);
358a7e14dcfSSatish Balay 
359ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
360a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
361a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
362a7e14dcfSSatish Balay           tao->ksp_its+=its;
3632d9aa51bSJason Sarich           tao->ksp_tot_its+=its;
364ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
365a7e14dcfSSatish Balay 
36653506e15SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
367a7e14dcfSSatish Balay         }
368a7e14dcfSSatish Balay       }
369a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
370a7e14dcfSSatish Balay       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
371a7e14dcfSSatish Balay       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&
372a7e14dcfSSatish Balay           (NTR_PC_BFGS == tr->pc_type) && (bfgsUpdates > 1)) {
373a7e14dcfSSatish Balay         /* Preconditioner is numerically indefinite; reset the
374a7e14dcfSSatish Balay            approximate if using BFGS preconditioning. */
375a7e14dcfSSatish Balay 
376*cd929ea3SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
377*cd929ea3SAlp Dener         ierr = MatLMVMSetJ0Scale(tr->M, delta);CHKERRQ(ierr);
378*cd929ea3SAlp Dener         ierr = MatLMVMReset(tr->M, PETSC_FALSE);CHKERRQ(ierr);
379a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
380a7e14dcfSSatish Balay         bfgsUpdates = 1;
381a7e14dcfSSatish Balay       }
382a7e14dcfSSatish Balay 
383a7e14dcfSSatish Balay       if (NTR_UPDATE_REDUCTION == tr->update_type) {
384a7e14dcfSSatish Balay         /* Get predicted reduction */
385ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
386a7e14dcfSSatish Balay         if (prered >= 0.0) {
387a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
388a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
389a7e14dcfSSatish Balay              be rejected! */
390a7e14dcfSSatish Balay           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
391a7e14dcfSSatish Balay         }
392a7e14dcfSSatish Balay         else {
393a7e14dcfSSatish Balay           /* Compute trial step and function value */
394a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
395a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
396a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
399a7e14dcfSSatish Balay             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
400a7e14dcfSSatish Balay           } else {
401a7e14dcfSSatish Balay             /* Compute and actual reduction */
402a7e14dcfSSatish Balay             actred = f - ftrial;
403a7e14dcfSSatish Balay             prered = -prered;
404a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
405a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
406a7e14dcfSSatish Balay               kappa = 1.0;
407a7e14dcfSSatish Balay             }
408a7e14dcfSSatish Balay             else {
409a7e14dcfSSatish Balay               kappa = actred / prered;
410a7e14dcfSSatish Balay             }
411a7e14dcfSSatish Balay 
412a7e14dcfSSatish Balay             /* Accept or reject the step and update radius */
413a7e14dcfSSatish Balay             if (kappa < tr->eta1) {
414a7e14dcfSSatish Balay               /* Reject the step */
415a7e14dcfSSatish Balay               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
416a7e14dcfSSatish Balay             }
417a7e14dcfSSatish Balay             else {
418a7e14dcfSSatish Balay               /* Accept the step */
419a7e14dcfSSatish Balay               if (kappa < tr->eta2) {
420a7e14dcfSSatish Balay                 /* Marginal bad step */
421a7e14dcfSSatish Balay                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
422a7e14dcfSSatish Balay               }
423a7e14dcfSSatish Balay               else if (kappa < tr->eta3) {
424a7e14dcfSSatish Balay                 /* Reasonable step */
425a7e14dcfSSatish Balay                 tao->trust = tr->alpha3 * tao->trust;
426a7e14dcfSSatish Balay               }
427a7e14dcfSSatish Balay               else if (kappa < tr->eta4) {
428a7e14dcfSSatish Balay                 /* Good step */
429a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
430a7e14dcfSSatish Balay               }
431a7e14dcfSSatish Balay               else {
432a7e14dcfSSatish Balay                 /* Very good step */
433a7e14dcfSSatish Balay                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
434a7e14dcfSSatish Balay               }
435a7e14dcfSSatish Balay               break;
436a7e14dcfSSatish Balay             }
437a7e14dcfSSatish Balay           }
438a7e14dcfSSatish Balay         }
439a7e14dcfSSatish Balay       }
440a7e14dcfSSatish Balay       else {
441a7e14dcfSSatish Balay         /* Get predicted reduction */
442ba7fe8fbSTodd Munson         ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
443a7e14dcfSSatish Balay         if (prered >= 0.0) {
444a7e14dcfSSatish Balay           /* The predicted reduction has the wrong sign.  This cannot
445a7e14dcfSSatish Balay              happen in infinite precision arithmetic.  Step should
446a7e14dcfSSatish Balay              be rejected! */
447a7e14dcfSSatish Balay           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
448a7e14dcfSSatish Balay         }
449a7e14dcfSSatish Balay         else {
450a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
451a7e14dcfSSatish Balay           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
452a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
453a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
454a7e14dcfSSatish Balay             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
455a7e14dcfSSatish Balay           }
456a7e14dcfSSatish Balay           else {
457a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
458a7e14dcfSSatish Balay             actred = f - ftrial;
459a7e14dcfSSatish Balay             prered = -prered;
460a7e14dcfSSatish Balay             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
461a7e14dcfSSatish Balay                 (PetscAbsScalar(prered) <= tr->epsilon)) {
462a7e14dcfSSatish Balay               kappa = 1.0;
463a7e14dcfSSatish Balay             }
464a7e14dcfSSatish Balay             else {
465a7e14dcfSSatish Balay               kappa = actred / prered;
466a7e14dcfSSatish Balay             }
467a7e14dcfSSatish Balay 
468a7e14dcfSSatish Balay             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
469a7e14dcfSSatish Balay             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
470a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
471a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
472a7e14dcfSSatish Balay 
473a7e14dcfSSatish Balay             if (kappa >= 1.0 - tr->mu1) {
474a7e14dcfSSatish Balay               /* Great agreement; accept step and update radius */
475a7e14dcfSSatish Balay               if (tau_max < 1.0) {
476a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
477a7e14dcfSSatish Balay               }
478a7e14dcfSSatish Balay               else if (tau_max > tr->gamma4) {
479a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
480a7e14dcfSSatish Balay               }
481a7e14dcfSSatish Balay               else {
482a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
483a7e14dcfSSatish Balay               }
484a7e14dcfSSatish Balay               break;
485a7e14dcfSSatish Balay             }
486a7e14dcfSSatish Balay             else if (kappa >= 1.0 - tr->mu2) {
487a7e14dcfSSatish Balay               /* Good agreement */
488a7e14dcfSSatish Balay 
489a7e14dcfSSatish Balay               if (tau_max < tr->gamma2) {
490a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
491a7e14dcfSSatish Balay               }
492a7e14dcfSSatish Balay               else if (tau_max > tr->gamma3) {
493a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
494a7e14dcfSSatish Balay               }
495a7e14dcfSSatish Balay               else if (tau_max < 1.0) {
496a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
497a7e14dcfSSatish Balay               }
498a7e14dcfSSatish Balay               else {
499a7e14dcfSSatish Balay                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
500a7e14dcfSSatish Balay               }
501a7e14dcfSSatish Balay               break;
502a7e14dcfSSatish Balay             }
503a7e14dcfSSatish Balay             else {
504a7e14dcfSSatish Balay               /* Not good agreement */
505a7e14dcfSSatish Balay               if (tau_min > 1.0) {
506a7e14dcfSSatish Balay                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
507a7e14dcfSSatish Balay               }
508a7e14dcfSSatish Balay               else if (tau_max < tr->gamma1) {
509a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
510a7e14dcfSSatish Balay               }
511a7e14dcfSSatish Balay               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
512a7e14dcfSSatish Balay                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
513a7e14dcfSSatish Balay               }
514a7e14dcfSSatish Balay               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
515a7e14dcfSSatish Balay                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
516a7e14dcfSSatish Balay                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
517a7e14dcfSSatish Balay               }
518a7e14dcfSSatish Balay               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
519a7e14dcfSSatish Balay                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
520a7e14dcfSSatish Balay                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
521a7e14dcfSSatish Balay               }
522a7e14dcfSSatish Balay               else {
523a7e14dcfSSatish Balay                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
524a7e14dcfSSatish Balay               }
525a7e14dcfSSatish Balay             }
526a7e14dcfSSatish Balay           }
527a7e14dcfSSatish Balay         }
528a7e14dcfSSatish Balay       }
529a7e14dcfSSatish Balay 
530a7e14dcfSSatish Balay       /* The step computed was not good and the radius was decreased.
531a7e14dcfSSatish Balay          Monitor the radius to terminate. */
5323ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
5333ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,tao->trust);CHKERRQ(ierr);
5343ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
535a7e14dcfSSatish Balay     }
536a7e14dcfSSatish Balay 
537a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
538a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tr->max_radius);
539a7e14dcfSSatish Balay 
5403ecd9318SAlp Dener     if (tao->reason == TAO_CONTINUE_ITERATING) {
541a7e14dcfSSatish Balay       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
542a7e14dcfSSatish Balay       f = ftrial;
543302440fdSBarry Smith       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr);
544a9603a14SPatrick Farrell       ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
54553506e15SBarry Smith       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
546a7e14dcfSSatish Balay       needH = 1;
5473ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
5483ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,tao->trust);CHKERRQ(ierr);
5493ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
550a7e14dcfSSatish Balay     }
551a7e14dcfSSatish Balay   }
552a7e14dcfSSatish Balay   PetscFunctionReturn(0);
553a7e14dcfSSatish Balay }
554a7e14dcfSSatish Balay 
555a7e14dcfSSatish Balay /*------------------------------------------------------------*/
556441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao)
557a7e14dcfSSatish Balay {
558a7e14dcfSSatish Balay   TAO_NTR *tr = (TAO_NTR *)tao->data;
559a7e14dcfSSatish Balay   PetscErrorCode ierr;
560a7e14dcfSSatish Balay 
561a7e14dcfSSatish Balay   PetscFunctionBegin;
562a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
563a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
564a7e14dcfSSatish Balay   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
565a7e14dcfSSatish Balay 
566a7e14dcfSSatish Balay   tr->Diag = 0;
567a7e14dcfSSatish Balay   tr->M = 0;
568a7e14dcfSSatish Balay   PetscFunctionReturn(0);
569a7e14dcfSSatish Balay }
570a7e14dcfSSatish Balay 
571a7e14dcfSSatish Balay /*------------------------------------------------------------*/
572441846f8SBarry Smith static PetscErrorCode TaoDestroy_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   if (tao->setupcalled) {
579a7e14dcfSSatish Balay     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
580a7e14dcfSSatish Balay   }
581a7e14dcfSSatish Balay   ierr = MatDestroy(&tr->M);CHKERRQ(ierr);
582a7e14dcfSSatish Balay   ierr = VecDestroy(&tr->Diag);CHKERRQ(ierr);
583a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
584a7e14dcfSSatish Balay   PetscFunctionReturn(0);
585a7e14dcfSSatish Balay }
586a7e14dcfSSatish Balay 
587a7e14dcfSSatish Balay /*------------------------------------------------------------*/
5884416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(PetscOptionItems *PetscOptionsObject,Tao tao)
589a7e14dcfSSatish Balay {
590a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
591a7e14dcfSSatish Balay   PetscErrorCode ierr;
592a7e14dcfSSatish Balay 
593a7e14dcfSSatish Balay   PetscFunctionBegin;
5941a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
59594ae4db5SBarry 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);
59694ae4db5SBarry 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);
59794ae4db5SBarry 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);
59894ae4db5SBarry 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);
59994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr);
60094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr);
60194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr);
60294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr);
60394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr);
60494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr);
60594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr);
60694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr);
60794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr);
60894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr);
60994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr);
61094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr);
61194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr);
61294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr);
61394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr);
61494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr);
61594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr);
61694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr);
61794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr);
61894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr);
61994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr);
62094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr);
62194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr);
62294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr);
62394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr);
62494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr);
625a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
626a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
627a7e14dcfSSatish Balay   PetscFunctionReturn(0);
628a7e14dcfSSatish Balay }
629a7e14dcfSSatish Balay 
630a7e14dcfSSatish Balay /*------------------------------------------------------------*/
631441846f8SBarry Smith static PetscErrorCode TaoView_NTR(Tao tao, PetscViewer viewer)
632a7e14dcfSSatish Balay {
633a7e14dcfSSatish Balay   TAO_NTR        *tr = (TAO_NTR *)tao->data;
634a7e14dcfSSatish Balay   PetscErrorCode ierr;
635a7e14dcfSSatish Balay   PetscInt       nrejects;
636a7e14dcfSSatish Balay   PetscBool      isascii;
63753506e15SBarry Smith 
638a7e14dcfSSatish Balay   PetscFunctionBegin;
639a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
640a7e14dcfSSatish Balay   if (isascii) {
641a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
642a7e14dcfSSatish Balay     if (NTR_PC_BFGS == tr->pc_type && tr->M) {
643*cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(tr->M, &nrejects);CHKERRQ(ierr);
644a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr);
645a7e14dcfSSatish Balay     }
646a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
647a7e14dcfSSatish Balay   }
648a7e14dcfSSatish Balay   PetscFunctionReturn(0);
649a7e14dcfSSatish Balay }
650a7e14dcfSSatish Balay 
651a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6521522df2eSJason Sarich /*MC
6531522df2eSJason Sarich   TAONTR - Newton's method with trust region for unconstrained minimization.
6541522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system.
6551522df2eSJason Sarich   NTR expects a KSP solver with a trust region radius.
6561522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
6571522df2eSJason Sarich 
6581522df2eSJason Sarich   Options Database Keys:
659fb90e4d1STodd Munson + -tao_ntr_pc_type - "none","ahess","bfgs","petsc"
6601522df2eSJason Sarich . -tao_ntr_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs"
6611522df2eSJason Sarich . -tao_ntr_init_type - "constant","direction","interpolation"
6621522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation"
6631522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius
6641522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius
6651522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
6661522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor
6671522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor
6681522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor
6691522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor
6701522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor
6711522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor
6721522df2eSJason Sarich . -tao_ntr_theta_i - thetha1 interpolation init factor
6731522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor
6741522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor
6751522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor
6761522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor
6771522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor
6781522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor
6791522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor
6801522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
6811522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
6821522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update
6831522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update
6841522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update
6851522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update
6861522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update
6871522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update
6881522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update
6891522df2eSJason Sarich 
6901eb8069cSJason Sarich   Level: beginner
6911522df2eSJason Sarich M*/
6921522df2eSJason Sarich 
693728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
694a7e14dcfSSatish Balay {
695a7e14dcfSSatish Balay   TAO_NTR *tr;
696a7e14dcfSSatish Balay   PetscErrorCode ierr;
697a7e14dcfSSatish Balay 
698a7e14dcfSSatish Balay   PetscFunctionBegin;
699a7e14dcfSSatish Balay 
7003c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
701a7e14dcfSSatish Balay 
702a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTR;
703a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTR;
704a7e14dcfSSatish Balay   tao->ops->view = TaoView_NTR;
705a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
706a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTR;
707a7e14dcfSSatish Balay 
7086552cf8aSJason Sarich   /* Override default settings (unless already changed) */
7096552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
7106552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
711a7e14dcfSSatish Balay   tao->data = (void*)tr;
712a7e14dcfSSatish Balay 
713a7e14dcfSSatish Balay   /*  Standard trust region update parameters */
714a7e14dcfSSatish Balay   tr->eta1 = 1.0e-4;
715a7e14dcfSSatish Balay   tr->eta2 = 0.25;
716a7e14dcfSSatish Balay   tr->eta3 = 0.50;
717a7e14dcfSSatish Balay   tr->eta4 = 0.90;
718a7e14dcfSSatish Balay 
719a7e14dcfSSatish Balay   tr->alpha1 = 0.25;
720a7e14dcfSSatish Balay   tr->alpha2 = 0.50;
721a7e14dcfSSatish Balay   tr->alpha3 = 1.00;
722a7e14dcfSSatish Balay   tr->alpha4 = 2.00;
723a7e14dcfSSatish Balay   tr->alpha5 = 4.00;
724a7e14dcfSSatish Balay 
725a7e14dcfSSatish Balay   /*  Interpolation trust region update parameters */
726a7e14dcfSSatish Balay   tr->mu1 = 0.10;
727a7e14dcfSSatish Balay   tr->mu2 = 0.50;
728a7e14dcfSSatish Balay 
729a7e14dcfSSatish Balay   tr->gamma1 = 0.25;
730a7e14dcfSSatish Balay   tr->gamma2 = 0.50;
731a7e14dcfSSatish Balay   tr->gamma3 = 2.00;
732a7e14dcfSSatish Balay   tr->gamma4 = 4.00;
733a7e14dcfSSatish Balay 
734a7e14dcfSSatish Balay   tr->theta = 0.05;
735a7e14dcfSSatish Balay 
736fb90e4d1STodd Munson   /*  Interpolation parameters for initialization */
737fb90e4d1STodd Munson   tr->mu1_i = 0.35;
738fb90e4d1STodd Munson   tr->mu2_i = 0.50;
739fb90e4d1STodd Munson 
740fb90e4d1STodd Munson   tr->gamma1_i = 0.0625;
741fb90e4d1STodd Munson   tr->gamma2_i = 0.50;
742fb90e4d1STodd Munson   tr->gamma3_i = 2.00;
743fb90e4d1STodd Munson   tr->gamma4_i = 5.00;
744fb90e4d1STodd Munson 
745fb90e4d1STodd Munson   tr->theta_i = 0.25;
746fb90e4d1STodd Munson 
747a7e14dcfSSatish Balay   tr->min_radius = 1.0e-10;
748a7e14dcfSSatish Balay   tr->max_radius = 1.0e10;
749a7e14dcfSSatish Balay   tr->epsilon    = 1.0e-6;
750a7e14dcfSSatish Balay 
751a7e14dcfSSatish Balay   tr->pc_type         = NTR_PC_BFGS;
752a7e14dcfSSatish Balay   tr->bfgs_scale_type = BFGS_SCALE_AHESS;
753a7e14dcfSSatish Balay   tr->init_type       = NTR_INIT_INTERPOLATION;
754a7e14dcfSSatish Balay   tr->update_type     = NTR_UPDATE_REDUCTION;
755a7e14dcfSSatish Balay 
756a7e14dcfSSatish Balay   /* Set linear solver to default for trust region */
757a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
75863b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
7595d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
760fb90e4d1STodd Munson   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
761a7e14dcfSSatish Balay   PetscFunctionReturn(0);
762a7e14dcfSSatish Balay }
763