xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision 63b1541575edd3f8a34e081896782a8d7a30de53)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
31daac38eSTodd Munson #include <../src/tao/unconstrained/impls/nls/nlsimpl.h>
4a7e14dcfSSatish Balay 
5aaa7dc30SBarry Smith #include <petscksp.h>
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay #define NLS_PC_NONE     0
8a7e14dcfSSatish Balay #define NLS_PC_AHESS    1
9a7e14dcfSSatish Balay #define NLS_PC_BFGS     2
10a7e14dcfSSatish Balay #define NLS_PC_PETSC    3
11a7e14dcfSSatish Balay #define NLS_PC_TYPES    4
12a7e14dcfSSatish Balay 
13a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS        0
14a7e14dcfSSatish Balay #define BFGS_SCALE_PHESS        1
15a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS         2
16a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES        3
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT         0
19a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION        1
20a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION    2
21a7e14dcfSSatish Balay #define NLS_INIT_TYPES            3
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay #define NLS_UPDATE_STEP           0
24a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION      1
25a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION  2
26a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES          3
27a7e14dcfSSatish Balay 
2887f595a5SBarry Smith static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"};
29a7e14dcfSSatish Balay 
3087f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
31a7e14dcfSSatish Balay 
3287f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"};
33a7e14dcfSSatish Balay 
3487f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"};
35a7e14dcfSSatish Balay 
361daac38eSTodd Munson /* Routine for BFGS preconditioner */
37a7e14dcfSSatish Balay 
381daac38eSTodd Munson static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
391daac38eSTodd Munson {
401daac38eSTodd Munson   PetscErrorCode ierr;
411daac38eSTodd Munson   Mat            M;
42a7e14dcfSSatish Balay 
431daac38eSTodd Munson   PetscFunctionBegin;
441daac38eSTodd Munson   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
451daac38eSTodd Munson   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
461daac38eSTodd Munson   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
471daac38eSTodd Munson   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
481daac38eSTodd Munson   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
491daac38eSTodd Munson   PetscFunctionReturn(0);
501daac38eSTodd Munson }
511daac38eSTodd Munson 
521daac38eSTodd Munson /*
53a7e14dcfSSatish Balay  Implements Newton's Method with a line search approach for solving
54a7e14dcfSSatish Balay  unconstrained minimization problems.  A More'-Thuente line search
55a7e14dcfSSatish Balay  is used to guarantee that the bfgs preconditioner remains positive
56a7e14dcfSSatish Balay  definite.
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay  The method can shift the Hessian matrix.  The shifting procedure is
59a7e14dcfSSatish Balay  adapted from the PATH algorithm for solving complementarity
60a7e14dcfSSatish Balay  problems.
61a7e14dcfSSatish Balay 
62a7e14dcfSSatish Balay  The linear system solve should be done with a conjugate gradient
631daac38eSTodd Munson  method, although any method can be used.
641daac38eSTodd Munson */
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay #define NLS_NEWTON              0
67a7e14dcfSSatish Balay #define NLS_BFGS                1
68a7e14dcfSSatish Balay #define NLS_SCALED_GRADIENT     2
69a7e14dcfSSatish Balay #define NLS_GRADIENT            3
70a7e14dcfSSatish Balay 
71441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao tao)
72a7e14dcfSSatish Balay {
73a7e14dcfSSatish Balay   PetscErrorCode               ierr;
74a7e14dcfSSatish Balay   TAO_NLS                      *nlsP = (TAO_NLS *)tao->data;
751daac38eSTodd Munson   KSPType                      ksp_type;
761daac38eSTodd Munson   PetscBool                    is_nash,is_stcg,is_gltr;
77a7e14dcfSSatish Balay   KSPConvergedReason           ksp_reason;
781daac38eSTodd Munson   PC                           pc;
79e4cb33bbSBarry Smith   TaoConvergedReason           reason;
801daac38eSTodd Munson   TaoLineSearchConvergedReason ls_reason;
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay   PetscReal                    fmin, ftrial, f_full, prered, actred, kappa, sigma;
83a7e14dcfSSatish Balay   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
84a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm, pert;
85a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
86a7e14dcfSSatish Balay   PetscReal                    delta;
87a7e14dcfSSatish Balay   PetscReal                    norm_d = 0.0, e_min;
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay   PetscInt                     stepType;
90a7e14dcfSSatish Balay   PetscInt                     bfgsUpdates = 0;
91a7e14dcfSSatish Balay   PetscInt                     n,N,kspits;
92b4690188SBarry Smith   PetscInt                     needH = 1;
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   PetscInt                     i_max = 5;
95a7e14dcfSSatish Balay   PetscInt                     j_max = 1;
96a7e14dcfSSatish Balay   PetscInt                     i, j;
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   PetscFunctionBegin;
99a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
100a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
101a7e14dcfSSatish Balay   }
102a7e14dcfSSatish Balay 
103a7e14dcfSSatish Balay   /* Initialized variables */
104a7e14dcfSSatish Balay   pert = nlsP->sval;
105a7e14dcfSSatish Balay 
1062d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
107a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
108a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
109a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
110a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
111a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
112a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
113a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
116ba7fe8fbSTodd Munson      Command automatically ignored for other methods
117ba7fe8fbSTodd Munson      Will be reset during the first iteration
118ba7fe8fbSTodd Munson   */
1191daac38eSTodd Munson   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
1201daac38eSTodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&is_nash);CHKERRQ(ierr);
1211daac38eSTodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&is_stcg);CHKERRQ(ierr);
1221daac38eSTodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&is_gltr);CHKERRQ(ierr);
1231daac38eSTodd Munson 
124ba7fe8fbSTodd Munson   ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
125a7e14dcfSSatish Balay 
1261daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
1271daac38eSTodd Munson     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
128a7e14dcfSSatish Balay     tao->trust = tao->trust0;
129a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
130a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
131a7e14dcfSSatish Balay   }
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay   /* Get vectors we will need */
134a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
135a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
136a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
137a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
138a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr);
139a7e14dcfSSatish Balay   }
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay   /* Check convergence criteria */
142a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
143a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
14487f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
145a7e14dcfSSatish Balay 
1468931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
14787f595a5SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay   /* create vectors for the limited memory preconditioner */
15087f595a5SBarry Smith   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
151a7e14dcfSSatish Balay     if (!nlsP->Diag) {
152a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
153a7e14dcfSSatish Balay     }
154a7e14dcfSSatish Balay   }
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay   /* Modify the preconditioner to use the bfgs approximation */
157a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
158a7e14dcfSSatish Balay   switch(nlsP->pc_type) {
159a7e14dcfSSatish Balay   case NLS_PC_NONE:
160a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1611a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
162a7e14dcfSSatish Balay     break;
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay   case NLS_PC_AHESS:
165a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1661a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
167baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
168a7e14dcfSSatish Balay     break;
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay   case NLS_PC_BFGS:
171a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
1721a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
173a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
174a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
175a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
176a7e14dcfSSatish Balay     break;
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay   default:
179a7e14dcfSSatish Balay     /* Use the pc method set by pc_type */
180a7e14dcfSSatish Balay     break;
181a7e14dcfSSatish Balay   }
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
184a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
1851daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
186a7e14dcfSSatish Balay     switch(nlsP->init_type) {
187a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
188a7e14dcfSSatish Balay       /* Use the initial radius specified */
189a7e14dcfSSatish Balay       break;
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
192a7e14dcfSSatish Balay       /* Use the initial radius specified */
193a7e14dcfSSatish Balay       max_radius = 0.0;
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
196a7e14dcfSSatish Balay         fmin = f;
197a7e14dcfSSatish Balay         sigma = 0.0;
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay         if (needH) {
200ffad9901SBarry Smith           ierr  = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
201a7e14dcfSSatish Balay           needH = 0;
202a7e14dcfSSatish Balay         }
203a7e14dcfSSatish Balay 
204a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
205a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
206a7e14dcfSSatish Balay           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
207a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
208a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
209a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
21087f595a5SBarry Smith           } else {
211a7e14dcfSSatish Balay             if (ftrial < fmin) {
212a7e14dcfSSatish Balay               fmin = ftrial;
213a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
214a7e14dcfSSatish Balay             }
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
217a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
220a7e14dcfSSatish Balay             actred = f - ftrial;
22187f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
222a7e14dcfSSatish Balay               kappa = 1.0;
22387f595a5SBarry Smith             } else {
224a7e14dcfSSatish Balay               kappa = actred / prered;
225a7e14dcfSSatish Balay             }
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
228a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
229a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
230a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
231a7e14dcfSSatish Balay 
232a7e14dcfSSatish Balay             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
233a7e14dcfSSatish Balay               /* Great agreement */
234a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
235a7e14dcfSSatish Balay 
236a7e14dcfSSatish Balay               if (tau_max < 1.0) {
237a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
23887f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
239a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
24087f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
241a7e14dcfSSatish Balay                 tau = tau_1;
24287f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
243a7e14dcfSSatish Balay                 tau = tau_2;
24487f595a5SBarry Smith               } else {
245a7e14dcfSSatish Balay                 tau = tau_max;
246a7e14dcfSSatish Balay               }
24787f595a5SBarry Smith             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
248a7e14dcfSSatish Balay               /* Good agreement */
249a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
250a7e14dcfSSatish Balay 
251a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
252a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
25387f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
254a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
25587f595a5SBarry Smith               } else {
256a7e14dcfSSatish Balay                 tau = tau_max;
257a7e14dcfSSatish Balay               }
25887f595a5SBarry Smith             } else {
259a7e14dcfSSatish Balay               /* Not good agreement */
260a7e14dcfSSatish Balay               if (tau_min > 1.0) {
261a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
26287f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
263a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
26487f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
265a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
26687f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
267a7e14dcfSSatish Balay                 tau = tau_1;
26887f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
269a7e14dcfSSatish Balay                 tau = tau_2;
27087f595a5SBarry Smith               } else {
271a7e14dcfSSatish Balay                 tau = tau_max;
272a7e14dcfSSatish Balay               }
273a7e14dcfSSatish Balay             }
274a7e14dcfSSatish Balay           }
275a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
276a7e14dcfSSatish Balay         }
277a7e14dcfSSatish Balay 
278a7e14dcfSSatish Balay         if (fmin < f) {
279a7e14dcfSSatish Balay           f = fmin;
280a7e14dcfSSatish Balay           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
281a7e14dcfSSatish Balay           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
282a7e14dcfSSatish Balay 
283a9603a14SPatrick Farrell           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
28487f595a5SBarry Smith           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
285a7e14dcfSSatish Balay           needH = 1;
286a7e14dcfSSatish Balay 
2878931d482SJason Sarich           ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
28887f595a5SBarry Smith           if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
289a7e14dcfSSatish Balay         }
290a7e14dcfSSatish Balay       }
291a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
292a7e14dcfSSatish Balay 
293a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
294a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
295a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
296a7e14dcfSSatish Balay       break;
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay     default:
299a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
300a7e14dcfSSatish Balay       tao->trust = 0.0;
301a7e14dcfSSatish Balay       break;
302a7e14dcfSSatish Balay     }
303a7e14dcfSSatish Balay   }
304a7e14dcfSSatish Balay 
305a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
306a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
307a7e14dcfSSatish Balay      since the function value may have decreased */
308a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type) {
309a7e14dcfSSatish Balay     if (f != 0.0) {
310a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
31187f595a5SBarry Smith     } else {
312a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
313a7e14dcfSSatish Balay     }
314a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
315a7e14dcfSSatish Balay   }
316a7e14dcfSSatish Balay 
317a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
318a7e14dcfSSatish Balay   nlsP->newt = 0;
319a7e14dcfSSatish Balay   nlsP->bfgs = 0;
320a7e14dcfSSatish Balay   nlsP->sgrad = 0;
321a7e14dcfSSatish Balay   nlsP->grad = 0;
322a7e14dcfSSatish Balay 
323a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
324a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
3258931d482SJason Sarich     ++tao->niter;
326ae93cb3cSJason Sarich     tao->ksp_its=0;
327a7e14dcfSSatish Balay 
328a7e14dcfSSatish Balay     /* Compute the Hessian */
329a7e14dcfSSatish Balay     if (needH) {
330ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
331a7e14dcfSSatish Balay     }
332a7e14dcfSSatish Balay 
33387f595a5SBarry Smith     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
334a7e14dcfSSatish Balay       /* Obtain diagonal for the bfgs preconditioner  */
335a7e14dcfSSatish Balay       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
336a7e14dcfSSatish Balay       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
337a7e14dcfSSatish Balay       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
338a7e14dcfSSatish Balay       ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
339a7e14dcfSSatish Balay     }
340a7e14dcfSSatish Balay 
341a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
342a7e14dcfSSatish Balay     if (pert > 0) {
343302440fdSBarry Smith       ierr = MatShift(tao->hessian, pert);CHKERRQ(ierr);
344a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
345a7e14dcfSSatish Balay         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
346a7e14dcfSSatish Balay       }
347a7e14dcfSSatish Balay     }
348a7e14dcfSSatish Balay 
349a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type) {
350a7e14dcfSSatish Balay       if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) {
351a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner  */
352a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
353a7e14dcfSSatish Balay         ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
354a7e14dcfSSatish Balay         ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
355a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
356a7e14dcfSSatish Balay       }
357a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
358a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
359a7e14dcfSSatish Balay       ++bfgsUpdates;
360a7e14dcfSSatish Balay     }
361a7e14dcfSSatish Balay 
362a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
36323ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
3641daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
365ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
366a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
367a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
368a7e14dcfSSatish Balay       tao->ksp_its+=kspits;
369ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
370ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
371a7e14dcfSSatish Balay 
372a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
373a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
374a7e14dcfSSatish Balay         if (norm_d > 0.0) {
375a7e14dcfSSatish Balay           tao->trust = norm_d;
376a7e14dcfSSatish Balay 
377a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
378a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
379a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
38087f595a5SBarry Smith         } else {
381a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
382a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
383a7e14dcfSSatish Balay           tao->trust = tao->trust0;
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
386a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
387a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
388a7e14dcfSSatish Balay 
389ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
390a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
391a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
392a7e14dcfSSatish Balay           tao->ksp_its+=kspits;
393ae93cb3cSJason Sarich           tao->ksp_tot_its+=kspits;
394ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
395a7e14dcfSSatish Balay 
39687f595a5SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
397a7e14dcfSSatish Balay         }
398a7e14dcfSSatish Balay       }
39987f595a5SBarry Smith     } else {
400a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
401a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
402a7e14dcfSSatish Balay       tao->ksp_its += kspits;
403ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
404a7e14dcfSSatish Balay     }
405a7e14dcfSSatish Balay     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
406a7e14dcfSSatish Balay     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
40787f595a5SBarry Smith     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
408a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
409a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
410a7e14dcfSSatish Balay 
411a7e14dcfSSatish Balay       if (f != 0.0) {
412a7e14dcfSSatish Balay         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
41387f595a5SBarry Smith       } else {
414a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
415a7e14dcfSSatish Balay       }
416a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
417a7e14dcfSSatish Balay       ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
418a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
419a7e14dcfSSatish Balay       bfgsUpdates = 1;
420a7e14dcfSSatish Balay     }
421a7e14dcfSSatish Balay 
422a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
423a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
42487f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
425a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
42687f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
427a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
42887f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
429a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
43087f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
431a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
43287f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
433a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
43487f595a5SBarry Smith     } else {
435a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
436a7e14dcfSSatish Balay     }
437a7e14dcfSSatish Balay 
438a7e14dcfSSatish Balay     /* Check for success (descent direction) */
439a7e14dcfSSatish Balay     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
440a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
441a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
442a7e14dcfSSatish Balay          Update the perturbation for next time */
443a7e14dcfSSatish Balay       if (pert <= 0.0) {
444a7e14dcfSSatish Balay         /* Initialize the perturbation */
445a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4461daac38eSTodd Munson         if (is_gltr) {
447ba7fe8fbSTodd Munson           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
448a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
449a7e14dcfSSatish Balay         }
45087f595a5SBarry Smith       } else {
451a7e14dcfSSatish Balay         /* Increase the perturbation */
452a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
453a7e14dcfSSatish Balay       }
454a7e14dcfSSatish Balay 
455a7e14dcfSSatish Balay       if (NLS_PC_BFGS != nlsP->pc_type) {
456a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
457a7e14dcfSSatish Balay            Must use gradient direction in this case */
458a7e14dcfSSatish Balay         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
459a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
460a7e14dcfSSatish Balay         ++nlsP->grad;
461a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
46287f595a5SBarry Smith       } else {
463a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
464a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
465a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
466a7e14dcfSSatish Balay 
467a7e14dcfSSatish Balay         /* Check for success (descent direction) */
468a7e14dcfSSatish Balay         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
469a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
470a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
471a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
472a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
473a7e14dcfSSatish Balay              which is guaranteed to be descent */
474a7e14dcfSSatish Balay 
475a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
476a7e14dcfSSatish Balay 
477a7e14dcfSSatish Balay           if (f != 0.0) {
478a7e14dcfSSatish Balay             delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
47987f595a5SBarry Smith           } else {
480a7e14dcfSSatish Balay             delta = 2.0 / (gnorm*gnorm);
481a7e14dcfSSatish Balay           }
482a7e14dcfSSatish Balay           ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
483a7e14dcfSSatish Balay           ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
484a7e14dcfSSatish Balay           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
485a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
486a7e14dcfSSatish Balay           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
487a7e14dcfSSatish Balay 
488a7e14dcfSSatish Balay           bfgsUpdates = 1;
489a7e14dcfSSatish Balay           ++nlsP->sgrad;
490a7e14dcfSSatish Balay           stepType = NLS_SCALED_GRADIENT;
49187f595a5SBarry Smith         } else {
492a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
493a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
494a7e14dcfSSatish Balay             ++nlsP->sgrad;
495a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
49687f595a5SBarry Smith           } else {
497a7e14dcfSSatish Balay             ++nlsP->bfgs;
498a7e14dcfSSatish Balay             stepType = NLS_BFGS;
499a7e14dcfSSatish Balay           }
500a7e14dcfSSatish Balay         }
501a7e14dcfSSatish Balay       }
50287f595a5SBarry Smith     } else {
503a7e14dcfSSatish Balay       /* Computed Newton step is descent */
504a7e14dcfSSatish Balay       switch (ksp_reason) {
505a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
506a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
507a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
508a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
509a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
510a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
511a7e14dcfSSatish Balay         if (pert <= 0.0) {
512a7e14dcfSSatish Balay           /* Initialize the perturbation */
513a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
5141daac38eSTodd Munson           if (is_gltr) {
515ba7fe8fbSTodd Munson             ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
516a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
517a7e14dcfSSatish Balay           }
51887f595a5SBarry Smith         } else {
519a7e14dcfSSatish Balay           /* Increase the perturbation */
520a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
521a7e14dcfSSatish Balay         }
522a7e14dcfSSatish Balay         break;
523a7e14dcfSSatish Balay 
524a7e14dcfSSatish Balay       default:
525a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
526a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
527a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
528a7e14dcfSSatish Balay           pert = 0.0;
529a7e14dcfSSatish Balay         }
530a7e14dcfSSatish Balay         break;
531a7e14dcfSSatish Balay       }
532a7e14dcfSSatish Balay 
533a7e14dcfSSatish Balay       ++nlsP->newt;
534a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
535a7e14dcfSSatish Balay     }
536a7e14dcfSSatish Balay 
537a7e14dcfSSatish Balay     /* Perform the linesearch */
538a7e14dcfSSatish Balay     fold = f;
539a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
540a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
541a7e14dcfSSatish Balay 
542a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
543a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
544a7e14dcfSSatish Balay 
54587f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
546a7e14dcfSSatish Balay       f = fold;
547a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
548a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
549a7e14dcfSSatish Balay 
550a7e14dcfSSatish Balay       switch(stepType) {
551a7e14dcfSSatish Balay       case NLS_NEWTON:
552a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
553a7e14dcfSSatish Balay            Update the perturbation for next time */
554a7e14dcfSSatish Balay         if (pert <= 0.0) {
555a7e14dcfSSatish Balay           /* Initialize the perturbation */
556a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
5571daac38eSTodd Munson           if (is_gltr) {
558ba7fe8fbSTodd Munson             ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
559a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
560a7e14dcfSSatish Balay           }
56187f595a5SBarry Smith         } else {
562a7e14dcfSSatish Balay           /* Increase the perturbation */
563a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
564a7e14dcfSSatish Balay         }
565a7e14dcfSSatish Balay 
566a7e14dcfSSatish Balay         if (NLS_PC_BFGS != nlsP->pc_type) {
567a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
568a7e14dcfSSatish Balay              Must use gradient direction in this case */
569a7e14dcfSSatish Balay           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
570a7e14dcfSSatish Balay           ++nlsP->grad;
571a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
57287f595a5SBarry Smith         } else {
573a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
574a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
575a7e14dcfSSatish Balay           /* Check for success (descent direction) */
576a7e14dcfSSatish Balay           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
577a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
578a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
579a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
580a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
581a7e14dcfSSatish Balay 
582a7e14dcfSSatish Balay             if (f != 0.0) {
583a7e14dcfSSatish Balay               delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
58487f595a5SBarry Smith             } else {
585a7e14dcfSSatish Balay               delta = 2.0 / (gnorm*gnorm);
586a7e14dcfSSatish Balay             }
587a7e14dcfSSatish Balay             ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
588a7e14dcfSSatish Balay             ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
589a7e14dcfSSatish Balay             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
590a7e14dcfSSatish Balay             ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
591a7e14dcfSSatish Balay 
592a7e14dcfSSatish Balay             bfgsUpdates = 1;
593a7e14dcfSSatish Balay             ++nlsP->sgrad;
594a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
59587f595a5SBarry Smith           } else {
596a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
597a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
598a7e14dcfSSatish Balay               ++nlsP->sgrad;
599a7e14dcfSSatish Balay               stepType = NLS_SCALED_GRADIENT;
60087f595a5SBarry Smith             } else {
601a7e14dcfSSatish Balay               ++nlsP->bfgs;
602a7e14dcfSSatish Balay               stepType = NLS_BFGS;
603a7e14dcfSSatish Balay             }
604a7e14dcfSSatish Balay           }
605a7e14dcfSSatish Balay         }
606a7e14dcfSSatish Balay         break;
607a7e14dcfSSatish Balay 
608a7e14dcfSSatish Balay       case NLS_BFGS:
609a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
610a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
611a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
612a7e14dcfSSatish Balay 
613a7e14dcfSSatish Balay         if (f != 0.0) {
614a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
61587f595a5SBarry Smith         } else {
616a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
617a7e14dcfSSatish Balay         }
618a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
619a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
620a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
621a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
622a7e14dcfSSatish Balay 
623a7e14dcfSSatish Balay         bfgsUpdates = 1;
624a7e14dcfSSatish Balay         ++nlsP->sgrad;
625a7e14dcfSSatish Balay         stepType = NLS_SCALED_GRADIENT;
626a7e14dcfSSatish Balay         break;
627a7e14dcfSSatish Balay 
628a7e14dcfSSatish Balay       case NLS_SCALED_GRADIENT:
629a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
630a7e14dcfSSatish Balay            The scaled gradient step did not produce a new iterate;
631a7e14dcfSSatish Balay            attemp to use the gradient direction.
632a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
633a7e14dcfSSatish Balay 
634a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr);
635a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr);
636a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
637a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
638a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
639a7e14dcfSSatish Balay 
640a7e14dcfSSatish Balay         bfgsUpdates = 1;
641a7e14dcfSSatish Balay         ++nlsP->grad;
642a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
643a7e14dcfSSatish Balay         break;
644a7e14dcfSSatish Balay       }
645a7e14dcfSSatish Balay       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
646a7e14dcfSSatish Balay 
647a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
648a7e14dcfSSatish Balay       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
649a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
650a7e14dcfSSatish Balay     }
651a7e14dcfSSatish Balay 
65287f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
653a7e14dcfSSatish Balay       /* Failed to find an improving point */
654a7e14dcfSSatish Balay       f = fold;
655a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
656a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
657a7e14dcfSSatish Balay       step = 0.0;
658a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
659a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
660a7e14dcfSSatish Balay       break;
661a7e14dcfSSatish Balay     }
662a7e14dcfSSatish Balay 
663a7e14dcfSSatish Balay     /* Update trust region radius */
6641daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
665a7e14dcfSSatish Balay       switch(nlsP->update_type) {
666a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
667a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
668a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
669a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
670a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
67187f595a5SBarry Smith           } else if (step < nlsP->nu2) {
672a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
673a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
67487f595a5SBarry Smith           } else if (step < nlsP->nu3) {
675a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
676a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
677a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
67887f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
679a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
680a7e14dcfSSatish Balay             }
68187f595a5SBarry Smith           } else if (step < nlsP->nu4) {
682a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
683a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
68487f595a5SBarry Smith           } else {
685a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
686a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
687a7e14dcfSSatish Balay           }
68887f595a5SBarry Smith         } else {
689a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
690a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
691a7e14dcfSSatish Balay         }
692a7e14dcfSSatish Balay         break;
693a7e14dcfSSatish Balay 
694a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
695a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
696a7e14dcfSSatish Balay           /*  Get predicted reduction */
697a7e14dcfSSatish Balay 
698ba7fe8fbSTodd Munson           ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
699a7e14dcfSSatish Balay           if (prered >= 0.0) {
700a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
701a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
702a7e14dcfSSatish Balay             /*  be rejected! */
703a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
70487f595a5SBarry Smith           } else {
705a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
706a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
70787f595a5SBarry Smith             } else {
708a7e14dcfSSatish Balay               /*  Compute and actual reduction */
709a7e14dcfSSatish Balay               actred = fold - f_full;
710a7e14dcfSSatish Balay               prered = -prered;
711a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
712a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
713a7e14dcfSSatish Balay                 kappa = 1.0;
71487f595a5SBarry Smith               } else {
715a7e14dcfSSatish Balay                 kappa = actred / prered;
716a7e14dcfSSatish Balay               }
717a7e14dcfSSatish Balay 
718a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
719a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
720a7e14dcfSSatish Balay                 /*  Very bad step */
721a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
72287f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
723a7e14dcfSSatish Balay                 /*  Marginal bad step */
724a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
72587f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
726a7e14dcfSSatish Balay                 /*  Reasonable step */
727a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
728a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
72987f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
730a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
731a7e14dcfSSatish Balay                 }
73287f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
733a7e14dcfSSatish Balay                 /*  Good step */
734a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
73587f595a5SBarry Smith               } else {
736a7e14dcfSSatish Balay                 /*  Very good step */
737a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
738a7e14dcfSSatish Balay               }
739a7e14dcfSSatish Balay             }
740a7e14dcfSSatish Balay           }
74187f595a5SBarry Smith         } else {
742a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
743a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
744a7e14dcfSSatish Balay         }
745a7e14dcfSSatish Balay         break;
746a7e14dcfSSatish Balay 
747a7e14dcfSSatish Balay       default:
748a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
749ba7fe8fbSTodd Munson           ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
750a7e14dcfSSatish Balay           if (prered >= 0.0) {
751a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
752a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
753a7e14dcfSSatish Balay             /*  be rejected! */
754a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
75587f595a5SBarry Smith           } else {
756a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
757a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
75887f595a5SBarry Smith             } else {
759a7e14dcfSSatish Balay               actred = fold - f_full;
760a7e14dcfSSatish Balay               prered = -prered;
76187f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
762a7e14dcfSSatish Balay                 kappa = 1.0;
76387f595a5SBarry Smith               } else {
764a7e14dcfSSatish Balay                 kappa = actred / prered;
765a7e14dcfSSatish Balay               }
766a7e14dcfSSatish Balay 
767a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
768a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
769a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
770a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
771a7e14dcfSSatish Balay 
772a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
773a7e14dcfSSatish Balay                 /*  Great agreement */
774a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
775a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
77687f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
777a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
77887f595a5SBarry Smith                 } else {
779a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
780a7e14dcfSSatish Balay                 }
78187f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
782a7e14dcfSSatish Balay                 /*  Good agreement */
783a7e14dcfSSatish Balay 
784a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
785a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
78687f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
787a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
78887f595a5SBarry Smith                 } else if (tau_max < 1.0) {
789a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
79087f595a5SBarry Smith                 } else {
791a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
792a7e14dcfSSatish Balay                 }
79387f595a5SBarry Smith               } else {
794a7e14dcfSSatish Balay                 /*  Not good agreement */
795a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
796a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
79787f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
798a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
79987f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
800a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
80187f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
802a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
80387f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
804a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
80587f595a5SBarry Smith                 } else {
806a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
807a7e14dcfSSatish Balay                 }
808a7e14dcfSSatish Balay               }
809a7e14dcfSSatish Balay             }
810a7e14dcfSSatish Balay           }
81187f595a5SBarry Smith         } else {
812a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
813a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
814a7e14dcfSSatish Balay         }
815a7e14dcfSSatish Balay         break;
816a7e14dcfSSatish Balay       }
817a7e14dcfSSatish Balay 
818a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
819a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
820a7e14dcfSSatish Balay     }
821a7e14dcfSSatish Balay 
822a7e14dcfSSatish Balay     /*  Check for termination */
823a9603a14SPatrick Farrell     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
82487f595a5SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
825a7e14dcfSSatish Balay     needH = 1;
8268931d482SJason Sarich     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
827a7e14dcfSSatish Balay   }
828a7e14dcfSSatish Balay   PetscFunctionReturn(0);
829a7e14dcfSSatish Balay }
830a7e14dcfSSatish Balay 
831a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
832441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao)
833a7e14dcfSSatish Balay {
834a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
835a7e14dcfSSatish Balay   PetscErrorCode ierr;
836a7e14dcfSSatish Balay 
837a7e14dcfSSatish Balay   PetscFunctionBegin;
838a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
839a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
840a7e14dcfSSatish Balay   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
841a7e14dcfSSatish Balay   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
842a7e14dcfSSatish Balay   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
843a7e14dcfSSatish Balay   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
844a7e14dcfSSatish Balay   nlsP->Diag = 0;
845a7e14dcfSSatish Balay   nlsP->M = 0;
846a7e14dcfSSatish Balay   PetscFunctionReturn(0);
847a7e14dcfSSatish Balay }
848a7e14dcfSSatish Balay 
849a7e14dcfSSatish Balay /*------------------------------------------------------------*/
850441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao)
851a7e14dcfSSatish Balay {
852a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
853a7e14dcfSSatish Balay   PetscErrorCode ierr;
854a7e14dcfSSatish Balay 
855a7e14dcfSSatish Balay   PetscFunctionBegin;
856a7e14dcfSSatish Balay   if (tao->setupcalled) {
857a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
858a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
859a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
860a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
861a7e14dcfSSatish Balay   }
862a7e14dcfSSatish Balay   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
863a7e14dcfSSatish Balay   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
864a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
865a7e14dcfSSatish Balay   PetscFunctionReturn(0);
866a7e14dcfSSatish Balay }
867a7e14dcfSSatish Balay 
868a7e14dcfSSatish Balay /*------------------------------------------------------------*/
8694416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptionItems *PetscOptionsObject,Tao tao)
870a7e14dcfSSatish Balay {
871a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
872a7e14dcfSSatish Balay   PetscErrorCode ierr;
873a7e14dcfSSatish Balay 
874a7e14dcfSSatish Balay   PetscFunctionBegin;
8751a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
876a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_nls_pc_type", "pc type", "", NLS_PC, NLS_PC_TYPES, NLS_PC[nlsP->pc_type], &nlsP->pc_type, 0);CHKERRQ(ierr);
877a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_nls_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[nlsP->bfgs_scale_type], &nlsP->bfgs_scale_type, 0);CHKERRQ(ierr);
878a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, 0);CHKERRQ(ierr);
879a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, 0);CHKERRQ(ierr);
88094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr);
88194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr);
88294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr);
88394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr);
88494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr);
88594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr);
88694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr);
88794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr);
88894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr);
88994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr);
89094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr);
89194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr);
89294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr);
89394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr);
89494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr);
89594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr);
89694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr);
89794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr);
89894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr);
89994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr);
90094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr);
90194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr);
90294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr);
90394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr);
90494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr);
90594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr);
90694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr);
90794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr);
90894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr);
90994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr);
91094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr);
91194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr);
91294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr);
91394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr);
91494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr);
91594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr);
91694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr);
91794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr);
91894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr);
91994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr);
92094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr);
92194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr);
92294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr);
92394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr);
92494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr);
925a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
926a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
927a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
928a7e14dcfSSatish Balay   PetscFunctionReturn(0);
929a7e14dcfSSatish Balay }
930a7e14dcfSSatish Balay 
931a7e14dcfSSatish Balay 
932a7e14dcfSSatish Balay /*------------------------------------------------------------*/
933441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
934a7e14dcfSSatish Balay {
935a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
936a7e14dcfSSatish Balay   PetscInt       nrejects;
937a7e14dcfSSatish Balay   PetscBool      isascii;
938a7e14dcfSSatish Balay   PetscErrorCode ierr;
939a7e14dcfSSatish Balay 
940a7e14dcfSSatish Balay   PetscFunctionBegin;
941a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
942a7e14dcfSSatish Balay   if (isascii) {
943a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
944a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
945a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr);
946a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
947a7e14dcfSSatish Balay     }
948a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
949a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
950a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
951a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
952a7e14dcfSSatish Balay 
953a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
954a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
955a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
956a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
957a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
958a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
959a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
960a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
961a7e14dcfSSatish Balay   }
962a7e14dcfSSatish Balay   PetscFunctionReturn(0);
963a7e14dcfSSatish Balay }
964a7e14dcfSSatish Balay 
965a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
9664aa34175SJason Sarich /*MC
9674aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
9684aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
9694aa34175SJason Sarich   system of equations to obtain the step diretion dk:
9704aa34175SJason Sarich               Hk dk = -gk
9714aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
9724aa34175SJason Sarich   solve
9734aa34175SJason Sarich               min_t f(xk + t d_k)
9744aa34175SJason Sarich 
9754aa34175SJason Sarich     Options Database Keys:
9761daac38eSTodd Munson + -tao_nls_pc_type - "none","ahess","bfgs","petsc"
9774aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs"
9784aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation"
9794aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
9804aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
9814aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
9824aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
9834aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
9844aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
9854aa34175SJason Sarich . -tao_nls_pgfac - growth factor
9864aa34175SJason Sarich . -tao_nls_psfac - shrink factor
9874aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
9884aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
9894aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
9904aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
9914aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
9924aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
9934aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
9944aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
9954aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
9964aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
9974aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
9984aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
9994aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
10004aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
10014aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
10024aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
10034aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
10044aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
10054aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
10064aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
10074aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
10084aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
10094aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
10104aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
10111522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
10121522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
10131522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
10141522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
10151522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
10161522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
10171522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
10181eb8069cSJason Sarich 
10191eb8069cSJason Sarich   Level: beginner
10201522df2eSJason Sarich M*/
10214aa34175SJason Sarich 
1022728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
1023a7e14dcfSSatish Balay {
1024a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
10258caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
1026a7e14dcfSSatish Balay   PetscErrorCode ierr;
1027a7e14dcfSSatish Balay 
1028a7e14dcfSSatish Balay   PetscFunctionBegin;
10293c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1030a7e14dcfSSatish Balay 
1031a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
1032a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
1033a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
1034a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1035a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
1036a7e14dcfSSatish Balay 
10376552cf8aSJason Sarich   /* Override default settings (unless already changed) */
10386552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
10396552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
10406552cf8aSJason Sarich 
1041a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
1042a7e14dcfSSatish Balay 
1043a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
1044a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
1045a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
1046a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
1047a7e14dcfSSatish Balay 
1048a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
1049a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
1050a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
1051a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
1052a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
1053a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
1054a7e14dcfSSatish Balay 
1055a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
1056a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
1057a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
1058a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
1059a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
1060a7e14dcfSSatish Balay 
1061a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
1062a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
1063a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
1064a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
1065a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
1066a7e14dcfSSatish Balay 
1067a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
1068a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
1069a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
1070a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
1071a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
1072a7e14dcfSSatish Balay 
1073a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
1074a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
1075a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
1076a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
1077a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
1078a7e14dcfSSatish Balay 
1079a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
1080a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
1081a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
1082a7e14dcfSSatish Balay 
1083a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
1084a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
1085a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
1086a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
1087a7e14dcfSSatish Balay 
1088a7e14dcfSSatish Balay   nlsP->theta = 0.05;
1089a7e14dcfSSatish Balay 
1090a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
1091a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
1092a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
1093a7e14dcfSSatish Balay 
1094a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
1095a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
1096a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
1097a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
1098a7e14dcfSSatish Balay 
1099a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
1100a7e14dcfSSatish Balay 
1101a7e14dcfSSatish Balay   /*  Remaining parameters */
1102a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
1103a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
1104a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
1105a7e14dcfSSatish Balay 
1106a7e14dcfSSatish Balay   nlsP->pc_type         = NLS_PC_BFGS;
1107a7e14dcfSSatish Balay   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1108a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1109a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
1110a7e14dcfSSatish Balay 
1111a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1112*63b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1113a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1114441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
11155d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1116a7e14dcfSSatish Balay 
1117a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
1118a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1119*63b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
11205d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
11211daac38eSTodd Munson   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1122a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1123a7e14dcfSSatish Balay }
1124