xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision cd929ea3f739fd9f7b6394f772cb40b9d4e6d97c)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
21daac38eSTodd Munson #include <../src/tao/unconstrained/impls/nls/nlsimpl.h>
3a7e14dcfSSatish Balay 
4aaa7dc30SBarry Smith #include <petscksp.h>
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay #define NLS_PC_NONE     0
7a7e14dcfSSatish Balay #define NLS_PC_AHESS    1
8a7e14dcfSSatish Balay #define NLS_PC_BFGS     2
9a7e14dcfSSatish Balay #define NLS_PC_PETSC    3
10a7e14dcfSSatish Balay #define NLS_PC_TYPES    4
11a7e14dcfSSatish Balay 
12a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS        0
13a7e14dcfSSatish Balay #define BFGS_SCALE_PHESS        1
14a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS         2
15a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES        3
16a7e14dcfSSatish Balay 
17a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT         0
18a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION        1
19a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION    2
20a7e14dcfSSatish Balay #define NLS_INIT_TYPES            3
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay #define NLS_UPDATE_STEP           0
23a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION      1
24a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION  2
25a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES          3
26a7e14dcfSSatish Balay 
2787f595a5SBarry Smith static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"};
28a7e14dcfSSatish Balay 
2987f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
30a7e14dcfSSatish Balay 
3187f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"};
32a7e14dcfSSatish Balay 
3387f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"};
34a7e14dcfSSatish Balay 
35*cd929ea3SAlp Dener PetscErrorCode TaoNLSPreconBFGS(PC BFGSpc, Vec X, Vec Y)
361daac38eSTodd Munson {
371daac38eSTodd Munson   PetscErrorCode ierr;
38*cd929ea3SAlp Dener   Mat *M;
39a7e14dcfSSatish Balay 
401daac38eSTodd Munson   PetscFunctionBegin;
41*cd929ea3SAlp Dener   ierr = PCShellGetContext(BFGSpc, (void**)&M);CHKERRQ(ierr);
42*cd929ea3SAlp Dener   ierr = MatSolve(*M, X, Y);CHKERRQ(ierr);
431daac38eSTodd Munson   PetscFunctionReturn(0);
441daac38eSTodd Munson }
451daac38eSTodd Munson 
461daac38eSTodd Munson /*
47a7e14dcfSSatish Balay  Implements Newton's Method with a line search approach for solving
48a7e14dcfSSatish Balay  unconstrained minimization problems.  A More'-Thuente line search
49a7e14dcfSSatish Balay  is used to guarantee that the bfgs preconditioner remains positive
50a7e14dcfSSatish Balay  definite.
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay  The method can shift the Hessian matrix.  The shifting procedure is
53a7e14dcfSSatish Balay  adapted from the PATH algorithm for solving complementarity
54a7e14dcfSSatish Balay  problems.
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay  The linear system solve should be done with a conjugate gradient
571daac38eSTodd Munson  method, although any method can be used.
581daac38eSTodd Munson */
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay #define NLS_NEWTON              0
61a7e14dcfSSatish Balay #define NLS_BFGS                1
62a7e14dcfSSatish Balay #define NLS_SCALED_GRADIENT     2
63a7e14dcfSSatish Balay #define NLS_GRADIENT            3
64a7e14dcfSSatish Balay 
65441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao tao)
66a7e14dcfSSatish Balay {
67a7e14dcfSSatish Balay   PetscErrorCode               ierr;
68a7e14dcfSSatish Balay   TAO_NLS                      *nlsP = (TAO_NLS *)tao->data;
691daac38eSTodd Munson   KSPType                      ksp_type;
701daac38eSTodd Munson   PetscBool                    is_nash,is_stcg,is_gltr;
71a7e14dcfSSatish Balay   KSPConvergedReason           ksp_reason;
721daac38eSTodd Munson   PC                           pc;
731daac38eSTodd Munson   TaoLineSearchConvergedReason ls_reason;
74a7e14dcfSSatish Balay 
75a7e14dcfSSatish Balay   PetscReal                    fmin, ftrial, f_full, prered, actred, kappa, sigma;
76a7e14dcfSSatish Balay   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
77a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm, pert;
78a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
79a7e14dcfSSatish Balay   PetscReal                    delta;
80a7e14dcfSSatish Balay   PetscReal                    norm_d = 0.0, e_min;
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay   PetscInt                     stepType;
83a7e14dcfSSatish Balay   PetscInt                     bfgsUpdates = 0;
84a7e14dcfSSatish Balay   PetscInt                     n,N,kspits;
85b4690188SBarry Smith   PetscInt                     needH = 1;
86a7e14dcfSSatish Balay 
87a7e14dcfSSatish Balay   PetscInt                     i_max = 5;
88a7e14dcfSSatish Balay   PetscInt                     j_max = 1;
89a7e14dcfSSatish Balay   PetscInt                     i, j;
90a7e14dcfSSatish Balay 
91a7e14dcfSSatish Balay   PetscFunctionBegin;
92a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
93a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
94a7e14dcfSSatish Balay   }
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay   /* Initialized variables */
97a7e14dcfSSatish Balay   pert = nlsP->sval;
98a7e14dcfSSatish Balay 
992d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
100a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
101a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
102a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
103a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
104a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
105a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
106a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
107a7e14dcfSSatish Balay 
108a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
109ba7fe8fbSTodd Munson      Command automatically ignored for other methods
110ba7fe8fbSTodd Munson      Will be reset during the first iteration
111ba7fe8fbSTodd Munson   */
1121daac38eSTodd Munson   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
1131daac38eSTodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&is_nash);CHKERRQ(ierr);
1141daac38eSTodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&is_stcg);CHKERRQ(ierr);
1151daac38eSTodd Munson   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&is_gltr);CHKERRQ(ierr);
1161daac38eSTodd Munson 
117ba7fe8fbSTodd Munson   ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
118a7e14dcfSSatish Balay 
1191daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
1201daac38eSTodd Munson     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
121a7e14dcfSSatish Balay     tao->trust = tao->trust0;
122a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
123a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
124a7e14dcfSSatish Balay   }
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay   /* Get vectors we will need */
127a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
128a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
129a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
130*cd929ea3SAlp Dener     ierr = MatCreateLBFGS(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
131*cd929ea3SAlp Dener     ierr = MatLMVMAllocate(nlsP->M,tao->solution,tao->gradient);CHKERRQ(ierr);
132a7e14dcfSSatish Balay   }
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay   /* Check convergence criteria */
135a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
136a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
13787f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
138a7e14dcfSSatish Balay 
1393ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1403ecd9318SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1413ecd9318SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr);
1423ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1433ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   /* create vectors for the limited memory preconditioner */
14687f595a5SBarry Smith   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
147a7e14dcfSSatish Balay     if (!nlsP->Diag) {
148a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
149a7e14dcfSSatish Balay     }
150a7e14dcfSSatish Balay   }
151a7e14dcfSSatish Balay 
152a7e14dcfSSatish Balay   /* Modify the preconditioner to use the bfgs approximation */
153a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   switch(nlsP->pc_type) {
155a7e14dcfSSatish Balay   case NLS_PC_NONE:
156a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1571a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
158a7e14dcfSSatish Balay     break;
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay   case NLS_PC_AHESS:
161a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1621a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
163baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
164a7e14dcfSSatish Balay     break;
165a7e14dcfSSatish Balay 
166a7e14dcfSSatish Balay   case NLS_PC_BFGS:
167a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
1681a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
169a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
170a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
171*cd929ea3SAlp Dener     ierr = PCShellSetApply(pc, TaoNLSPreconBFGS);CHKERRQ(ierr);
172a7e14dcfSSatish Balay     break;
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay   default:
175a7e14dcfSSatish Balay     /* Use the pc method set by pc_type */
176a7e14dcfSSatish Balay     break;
177a7e14dcfSSatish Balay   }
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
180a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
1811daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
182a7e14dcfSSatish Balay     switch(nlsP->init_type) {
183a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
184a7e14dcfSSatish Balay       /* Use the initial radius specified */
185a7e14dcfSSatish Balay       break;
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
188a7e14dcfSSatish Balay       /* Use the initial radius specified */
189a7e14dcfSSatish Balay       max_radius = 0.0;
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
192a7e14dcfSSatish Balay         fmin = f;
193a7e14dcfSSatish Balay         sigma = 0.0;
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay         if (needH) {
196ffad9901SBarry Smith           ierr  = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
197a7e14dcfSSatish Balay           needH = 0;
198a7e14dcfSSatish Balay         }
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
201a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
202a7e14dcfSSatish Balay           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
203a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
204a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
205a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
20687f595a5SBarry Smith           } else {
207a7e14dcfSSatish Balay             if (ftrial < fmin) {
208a7e14dcfSSatish Balay               fmin = ftrial;
209a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
210a7e14dcfSSatish Balay             }
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
213a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
214a7e14dcfSSatish Balay 
215a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
216a7e14dcfSSatish Balay             actred = f - ftrial;
21787f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
218a7e14dcfSSatish Balay               kappa = 1.0;
21987f595a5SBarry Smith             } else {
220a7e14dcfSSatish Balay               kappa = actred / prered;
221a7e14dcfSSatish Balay             }
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
224a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
225a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
226a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
227a7e14dcfSSatish Balay 
228a7e14dcfSSatish Balay             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
229a7e14dcfSSatish Balay               /* Great agreement */
230a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
231a7e14dcfSSatish Balay 
232a7e14dcfSSatish Balay               if (tau_max < 1.0) {
233a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
23487f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
235a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
23687f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
237a7e14dcfSSatish Balay                 tau = tau_1;
23887f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
239a7e14dcfSSatish Balay                 tau = tau_2;
24087f595a5SBarry Smith               } else {
241a7e14dcfSSatish Balay                 tau = tau_max;
242a7e14dcfSSatish Balay               }
24387f595a5SBarry Smith             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
244a7e14dcfSSatish Balay               /* Good agreement */
245a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
248a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
24987f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
250a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
25187f595a5SBarry Smith               } else {
252a7e14dcfSSatish Balay                 tau = tau_max;
253a7e14dcfSSatish Balay               }
25487f595a5SBarry Smith             } else {
255a7e14dcfSSatish Balay               /* Not good agreement */
256a7e14dcfSSatish Balay               if (tau_min > 1.0) {
257a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
25887f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
259a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
26087f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
261a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
26287f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
263a7e14dcfSSatish Balay                 tau = tau_1;
26487f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
265a7e14dcfSSatish Balay                 tau = tau_2;
26687f595a5SBarry Smith               } else {
267a7e14dcfSSatish Balay                 tau = tau_max;
268a7e14dcfSSatish Balay               }
269a7e14dcfSSatish Balay             }
270a7e14dcfSSatish Balay           }
271a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
272a7e14dcfSSatish Balay         }
273a7e14dcfSSatish Balay 
274a7e14dcfSSatish Balay         if (fmin < f) {
275a7e14dcfSSatish Balay           f = fmin;
276a7e14dcfSSatish Balay           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
277a7e14dcfSSatish Balay           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
278a7e14dcfSSatish Balay 
279a9603a14SPatrick Farrell           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
28087f595a5SBarry Smith           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
281a7e14dcfSSatish Balay           needH = 1;
282a7e14dcfSSatish Balay 
2833ecd9318SAlp Dener           ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
2843ecd9318SAlp Dener           ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr);
2853ecd9318SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2863ecd9318SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
287a7e14dcfSSatish Balay         }
288a7e14dcfSSatish Balay       }
289a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
290a7e14dcfSSatish Balay 
291a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
292a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
293a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
294a7e14dcfSSatish Balay       break;
295a7e14dcfSSatish Balay 
296a7e14dcfSSatish Balay     default:
297a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
298a7e14dcfSSatish Balay       tao->trust = 0.0;
299a7e14dcfSSatish Balay       break;
300a7e14dcfSSatish Balay     }
301a7e14dcfSSatish Balay   }
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
304a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
305a7e14dcfSSatish Balay      since the function value may have decreased */
306a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type) {
307*cd929ea3SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
308*cd929ea3SAlp Dener     ierr = MatLMVMSetJ0Scale(nlsP->M, delta);CHKERRQ(ierr);
309a7e14dcfSSatish Balay   }
310a7e14dcfSSatish Balay 
311a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
312a7e14dcfSSatish Balay   nlsP->newt = 0;
313a7e14dcfSSatish Balay   nlsP->bfgs = 0;
314a7e14dcfSSatish Balay   nlsP->sgrad = 0;
315a7e14dcfSSatish Balay   nlsP->grad = 0;
316a7e14dcfSSatish Balay 
317a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
3183ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
3198931d482SJason Sarich     ++tao->niter;
320ae93cb3cSJason Sarich     tao->ksp_its=0;
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay     /* Compute the Hessian */
323a7e14dcfSSatish Balay     if (needH) {
324ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
325a7e14dcfSSatish Balay     }
326a7e14dcfSSatish Balay 
32787f595a5SBarry Smith     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
328a7e14dcfSSatish Balay       /* Obtain diagonal for the bfgs preconditioner  */
329a7e14dcfSSatish Balay       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
330a7e14dcfSSatish Balay       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
331a7e14dcfSSatish Balay       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
332*cd929ea3SAlp Dener       ierr = MatLMVMSetJ0Diag(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
333a7e14dcfSSatish Balay     }
334a7e14dcfSSatish Balay 
335a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
336a7e14dcfSSatish Balay     if (pert > 0) {
337302440fdSBarry Smith       ierr = MatShift(tao->hessian, pert);CHKERRQ(ierr);
338a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
339a7e14dcfSSatish Balay         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
340a7e14dcfSSatish Balay       }
341a7e14dcfSSatish Balay     }
342a7e14dcfSSatish Balay 
343a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type) {
344a7e14dcfSSatish Balay       if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) {
345a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner  */
346a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
347a7e14dcfSSatish Balay         ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
348a7e14dcfSSatish Balay         ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
349*cd929ea3SAlp Dener         ierr = MatLMVMSetJ0Diag(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
350a7e14dcfSSatish Balay       }
351a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
352a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
353a7e14dcfSSatish Balay       ++bfgsUpdates;
354a7e14dcfSSatish Balay     }
355a7e14dcfSSatish Balay 
356a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
35723ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
3581daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
359ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
360a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
361a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
362a7e14dcfSSatish Balay       tao->ksp_its+=kspits;
363ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
364ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
365a7e14dcfSSatish Balay 
366a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
367a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
368a7e14dcfSSatish Balay         if (norm_d > 0.0) {
369a7e14dcfSSatish Balay           tao->trust = norm_d;
370a7e14dcfSSatish Balay 
371a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
372a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
373a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
37487f595a5SBarry Smith         } else {
375a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
376a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
377a7e14dcfSSatish Balay           tao->trust = tao->trust0;
378a7e14dcfSSatish Balay 
379a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
380a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
381a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
382a7e14dcfSSatish Balay 
383ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
384a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
385a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
386a7e14dcfSSatish Balay           tao->ksp_its+=kspits;
387ae93cb3cSJason Sarich           tao->ksp_tot_its+=kspits;
388ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
389a7e14dcfSSatish Balay 
39087f595a5SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
391a7e14dcfSSatish Balay         }
392a7e14dcfSSatish Balay       }
39387f595a5SBarry Smith     } else {
394a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
395a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
396a7e14dcfSSatish Balay       tao->ksp_its += kspits;
397ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
398a7e14dcfSSatish Balay     }
399a7e14dcfSSatish Balay     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
400a7e14dcfSSatish Balay     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
40187f595a5SBarry Smith     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
402a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
403a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
404a7e14dcfSSatish Balay 
405*cd929ea3SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
406*cd929ea3SAlp Dener       ierr = MatLMVMSetJ0Scale(nlsP->M, delta);CHKERRQ(ierr);
407*cd929ea3SAlp Dener       ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr);
408a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
409a7e14dcfSSatish Balay       bfgsUpdates = 1;
410a7e14dcfSSatish Balay     }
411a7e14dcfSSatish Balay 
412a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
413a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
41487f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
415a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
41687f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
417a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
41887f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
419a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
42087f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
421a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
42287f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
423a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
42487f595a5SBarry Smith     } else {
425a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
426a7e14dcfSSatish Balay     }
427a7e14dcfSSatish Balay 
428a7e14dcfSSatish Balay     /* Check for success (descent direction) */
429a7e14dcfSSatish Balay     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
430a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
431a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
432a7e14dcfSSatish Balay          Update the perturbation for next time */
433a7e14dcfSSatish Balay       if (pert <= 0.0) {
434a7e14dcfSSatish Balay         /* Initialize the perturbation */
435a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4361daac38eSTodd Munson         if (is_gltr) {
437ba7fe8fbSTodd Munson           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
438a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
439a7e14dcfSSatish Balay         }
44087f595a5SBarry Smith       } else {
441a7e14dcfSSatish Balay         /* Increase the perturbation */
442a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
443a7e14dcfSSatish Balay       }
444a7e14dcfSSatish Balay 
445a7e14dcfSSatish Balay       if (NLS_PC_BFGS != nlsP->pc_type) {
446a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
447a7e14dcfSSatish Balay            Must use gradient direction in this case */
448a7e14dcfSSatish Balay         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
449a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
450a7e14dcfSSatish Balay         ++nlsP->grad;
451a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
45287f595a5SBarry Smith       } else {
453a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
454*cd929ea3SAlp Dener         ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
455a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
456a7e14dcfSSatish Balay 
457a7e14dcfSSatish Balay         /* Check for success (descent direction) */
458a7e14dcfSSatish Balay         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
459a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
460a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
461a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
462a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
463a7e14dcfSSatish Balay              which is guaranteed to be descent */
464a7e14dcfSSatish Balay 
465a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
466a7e14dcfSSatish Balay 
467*cd929ea3SAlp Dener           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
468*cd929ea3SAlp Dener           ierr = MatLMVMSetJ0Scale(nlsP->M, delta);CHKERRQ(ierr);
469*cd929ea3SAlp Dener           ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr);
470a7e14dcfSSatish Balay           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
471*cd929ea3SAlp Dener           ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
472a7e14dcfSSatish Balay           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay           bfgsUpdates = 1;
475a7e14dcfSSatish Balay           ++nlsP->sgrad;
476a7e14dcfSSatish Balay           stepType = NLS_SCALED_GRADIENT;
47787f595a5SBarry Smith         } else {
478a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
479a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
480a7e14dcfSSatish Balay             ++nlsP->sgrad;
481a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
48287f595a5SBarry Smith           } else {
483a7e14dcfSSatish Balay             ++nlsP->bfgs;
484a7e14dcfSSatish Balay             stepType = NLS_BFGS;
485a7e14dcfSSatish Balay           }
486a7e14dcfSSatish Balay         }
487a7e14dcfSSatish Balay       }
48887f595a5SBarry Smith     } else {
489a7e14dcfSSatish Balay       /* Computed Newton step is descent */
490a7e14dcfSSatish Balay       switch (ksp_reason) {
491a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
492a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
493a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
494a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
495a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
496a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
497a7e14dcfSSatish Balay         if (pert <= 0.0) {
498a7e14dcfSSatish Balay           /* Initialize the perturbation */
499a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
5001daac38eSTodd Munson           if (is_gltr) {
501ba7fe8fbSTodd Munson             ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
502a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
503a7e14dcfSSatish Balay           }
50487f595a5SBarry Smith         } else {
505a7e14dcfSSatish Balay           /* Increase the perturbation */
506a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
507a7e14dcfSSatish Balay         }
508a7e14dcfSSatish Balay         break;
509a7e14dcfSSatish Balay 
510a7e14dcfSSatish Balay       default:
511a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
512a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
513a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
514a7e14dcfSSatish Balay           pert = 0.0;
515a7e14dcfSSatish Balay         }
516a7e14dcfSSatish Balay         break;
517a7e14dcfSSatish Balay       }
518a7e14dcfSSatish Balay 
519a7e14dcfSSatish Balay       ++nlsP->newt;
520a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
521a7e14dcfSSatish Balay     }
522a7e14dcfSSatish Balay 
523a7e14dcfSSatish Balay     /* Perform the linesearch */
524a7e14dcfSSatish Balay     fold = f;
525a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
526a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
527a7e14dcfSSatish Balay 
528a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
529a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
530a7e14dcfSSatish Balay 
53187f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
532a7e14dcfSSatish Balay       f = fold;
533a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
534a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
535a7e14dcfSSatish Balay 
536a7e14dcfSSatish Balay       switch(stepType) {
537a7e14dcfSSatish Balay       case NLS_NEWTON:
538a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
539a7e14dcfSSatish Balay            Update the perturbation for next time */
540a7e14dcfSSatish Balay         if (pert <= 0.0) {
541a7e14dcfSSatish Balay           /* Initialize the perturbation */
542a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
5431daac38eSTodd Munson           if (is_gltr) {
544ba7fe8fbSTodd Munson             ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
545a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
546a7e14dcfSSatish Balay           }
54787f595a5SBarry Smith         } else {
548a7e14dcfSSatish Balay           /* Increase the perturbation */
549a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
550a7e14dcfSSatish Balay         }
551a7e14dcfSSatish Balay 
552a7e14dcfSSatish Balay         if (NLS_PC_BFGS != nlsP->pc_type) {
553a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
554a7e14dcfSSatish Balay              Must use gradient direction in this case */
555a7e14dcfSSatish Balay           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
556a7e14dcfSSatish Balay           ++nlsP->grad;
557a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
55887f595a5SBarry Smith         } else {
559a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
560*cd929ea3SAlp Dener           ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
561a7e14dcfSSatish Balay           /* Check for success (descent direction) */
562a7e14dcfSSatish Balay           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
563a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
564a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
565a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
566a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
567a7e14dcfSSatish Balay 
568*cd929ea3SAlp Dener             delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
569*cd929ea3SAlp Dener             ierr = MatLMVMSetJ0Scale(nlsP->M, delta);CHKERRQ(ierr);
570*cd929ea3SAlp Dener             ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr);
571a7e14dcfSSatish Balay             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
572*cd929ea3SAlp Dener             ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
573a7e14dcfSSatish Balay 
574a7e14dcfSSatish Balay             bfgsUpdates = 1;
575a7e14dcfSSatish Balay             ++nlsP->sgrad;
576a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
57787f595a5SBarry Smith           } else {
578a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
579a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
580a7e14dcfSSatish Balay               ++nlsP->sgrad;
581a7e14dcfSSatish Balay               stepType = NLS_SCALED_GRADIENT;
58287f595a5SBarry Smith             } else {
583a7e14dcfSSatish Balay               ++nlsP->bfgs;
584a7e14dcfSSatish Balay               stepType = NLS_BFGS;
585a7e14dcfSSatish Balay             }
586a7e14dcfSSatish Balay           }
587a7e14dcfSSatish Balay         }
588a7e14dcfSSatish Balay         break;
589a7e14dcfSSatish Balay 
590a7e14dcfSSatish Balay       case NLS_BFGS:
591a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
592a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
593a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
594a7e14dcfSSatish Balay 
595*cd929ea3SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
596*cd929ea3SAlp Dener         ierr = MatLMVMSetJ0Scale(nlsP->M, delta);CHKERRQ(ierr);
597*cd929ea3SAlp Dener         ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr);
598a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
599*cd929ea3SAlp Dener         ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay         bfgsUpdates = 1;
602a7e14dcfSSatish Balay         ++nlsP->sgrad;
603a7e14dcfSSatish Balay         stepType = NLS_SCALED_GRADIENT;
604a7e14dcfSSatish Balay         break;
605a7e14dcfSSatish Balay 
606a7e14dcfSSatish Balay       case NLS_SCALED_GRADIENT:
607a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
608a7e14dcfSSatish Balay            The scaled gradient step did not produce a new iterate;
609a7e14dcfSSatish Balay            attemp to use the gradient direction.
610a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
611a7e14dcfSSatish Balay 
612*cd929ea3SAlp Dener         ierr = MatLMVMSetJ0Scale(nlsP->M,0);CHKERRQ(ierr);
613*cd929ea3SAlp Dener         ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr);
614a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
615*cd929ea3SAlp Dener         ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
616a7e14dcfSSatish Balay 
617a7e14dcfSSatish Balay         bfgsUpdates = 1;
618a7e14dcfSSatish Balay         ++nlsP->grad;
619a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
620a7e14dcfSSatish Balay         break;
621a7e14dcfSSatish Balay       }
622a7e14dcfSSatish Balay       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
623a7e14dcfSSatish Balay 
624a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
625a7e14dcfSSatish Balay       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
626a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
627a7e14dcfSSatish Balay     }
628a7e14dcfSSatish Balay 
62987f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
630a7e14dcfSSatish Balay       /* Failed to find an improving point */
631a7e14dcfSSatish Balay       f = fold;
632a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
633a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
634a7e14dcfSSatish Balay       step = 0.0;
635a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
636a7e14dcfSSatish Balay       break;
637a7e14dcfSSatish Balay     }
638a7e14dcfSSatish Balay 
639a7e14dcfSSatish Balay     /* Update trust region radius */
6401daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
641a7e14dcfSSatish Balay       switch(nlsP->update_type) {
642a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
643a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
644a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
645a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
646a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
64787f595a5SBarry Smith           } else if (step < nlsP->nu2) {
648a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
649a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
65087f595a5SBarry Smith           } else if (step < nlsP->nu3) {
651a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
652a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
653a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
65487f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
655a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
656a7e14dcfSSatish Balay             }
65787f595a5SBarry Smith           } else if (step < nlsP->nu4) {
658a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
659a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
66087f595a5SBarry Smith           } else {
661a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
662a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
663a7e14dcfSSatish Balay           }
66487f595a5SBarry Smith         } else {
665a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
666a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
667a7e14dcfSSatish Balay         }
668a7e14dcfSSatish Balay         break;
669a7e14dcfSSatish Balay 
670a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
671a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
672a7e14dcfSSatish Balay           /*  Get predicted reduction */
673a7e14dcfSSatish Balay 
674ba7fe8fbSTodd Munson           ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
675a7e14dcfSSatish Balay           if (prered >= 0.0) {
676a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
677a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
678a7e14dcfSSatish Balay             /*  be rejected! */
679a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
68087f595a5SBarry Smith           } else {
681a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
682a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
68387f595a5SBarry Smith             } else {
684a7e14dcfSSatish Balay               /*  Compute and actual reduction */
685a7e14dcfSSatish Balay               actred = fold - f_full;
686a7e14dcfSSatish Balay               prered = -prered;
687a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
688a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
689a7e14dcfSSatish Balay                 kappa = 1.0;
69087f595a5SBarry Smith               } else {
691a7e14dcfSSatish Balay                 kappa = actred / prered;
692a7e14dcfSSatish Balay               }
693a7e14dcfSSatish Balay 
694a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
695a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
696a7e14dcfSSatish Balay                 /*  Very bad step */
697a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
69887f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
699a7e14dcfSSatish Balay                 /*  Marginal bad step */
700a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
70187f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
702a7e14dcfSSatish Balay                 /*  Reasonable step */
703a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
704a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
70587f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
706a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
707a7e14dcfSSatish Balay                 }
70887f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
709a7e14dcfSSatish Balay                 /*  Good step */
710a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
71187f595a5SBarry Smith               } else {
712a7e14dcfSSatish Balay                 /*  Very good step */
713a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
714a7e14dcfSSatish Balay               }
715a7e14dcfSSatish Balay             }
716a7e14dcfSSatish Balay           }
71787f595a5SBarry Smith         } else {
718a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
719a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
720a7e14dcfSSatish Balay         }
721a7e14dcfSSatish Balay         break;
722a7e14dcfSSatish Balay 
723a7e14dcfSSatish Balay       default:
724a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
725ba7fe8fbSTodd Munson           ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
726a7e14dcfSSatish Balay           if (prered >= 0.0) {
727a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
728a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
729a7e14dcfSSatish Balay             /*  be rejected! */
730a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
73187f595a5SBarry Smith           } else {
732a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
733a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
73487f595a5SBarry Smith             } else {
735a7e14dcfSSatish Balay               actred = fold - f_full;
736a7e14dcfSSatish Balay               prered = -prered;
73787f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
738a7e14dcfSSatish Balay                 kappa = 1.0;
73987f595a5SBarry Smith               } else {
740a7e14dcfSSatish Balay                 kappa = actred / prered;
741a7e14dcfSSatish Balay               }
742a7e14dcfSSatish Balay 
743a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
744a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
745a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
746a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
747a7e14dcfSSatish Balay 
748a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
749a7e14dcfSSatish Balay                 /*  Great agreement */
750a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
751a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
75287f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
753a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
75487f595a5SBarry Smith                 } else {
755a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
756a7e14dcfSSatish Balay                 }
75787f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
758a7e14dcfSSatish Balay                 /*  Good agreement */
759a7e14dcfSSatish Balay 
760a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
761a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
76287f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
763a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
76487f595a5SBarry Smith                 } else if (tau_max < 1.0) {
765a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
76687f595a5SBarry Smith                 } else {
767a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
768a7e14dcfSSatish Balay                 }
76987f595a5SBarry Smith               } else {
770a7e14dcfSSatish Balay                 /*  Not good agreement */
771a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
772a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
77387f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
774a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
77587f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
776a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
77787f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
778a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
77987f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
780a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
78187f595a5SBarry Smith                 } else {
782a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
783a7e14dcfSSatish Balay                 }
784a7e14dcfSSatish Balay               }
785a7e14dcfSSatish Balay             }
786a7e14dcfSSatish Balay           }
78787f595a5SBarry Smith         } else {
788a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
789a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
790a7e14dcfSSatish Balay         }
791a7e14dcfSSatish Balay         break;
792a7e14dcfSSatish Balay       }
793a7e14dcfSSatish Balay 
794a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
795a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
796a7e14dcfSSatish Balay     }
797a7e14dcfSSatish Balay 
798a7e14dcfSSatish Balay     /*  Check for termination */
799a9603a14SPatrick Farrell     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
80087f595a5SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
801a7e14dcfSSatish Balay     needH = 1;
8023ecd9318SAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
8033ecd9318SAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr);
8043ecd9318SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
805a7e14dcfSSatish Balay   }
806a7e14dcfSSatish Balay   PetscFunctionReturn(0);
807a7e14dcfSSatish Balay }
808a7e14dcfSSatish Balay 
809a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
810441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao)
811a7e14dcfSSatish Balay {
812a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
813a7e14dcfSSatish Balay   PetscErrorCode ierr;
814a7e14dcfSSatish Balay 
815a7e14dcfSSatish Balay   PetscFunctionBegin;
816a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
817a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
818a7e14dcfSSatish Balay   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
819a7e14dcfSSatish Balay   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
820a7e14dcfSSatish Balay   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
821a7e14dcfSSatish Balay   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
822a7e14dcfSSatish Balay   nlsP->Diag = 0;
823a7e14dcfSSatish Balay   nlsP->M = 0;
824a7e14dcfSSatish Balay   PetscFunctionReturn(0);
825a7e14dcfSSatish Balay }
826a7e14dcfSSatish Balay 
827a7e14dcfSSatish Balay /*------------------------------------------------------------*/
828441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao)
829a7e14dcfSSatish Balay {
830a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
831a7e14dcfSSatish Balay   PetscErrorCode ierr;
832a7e14dcfSSatish Balay 
833a7e14dcfSSatish Balay   PetscFunctionBegin;
834a7e14dcfSSatish Balay   if (tao->setupcalled) {
835a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
836a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
837a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
838a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
839a7e14dcfSSatish Balay   }
840a7e14dcfSSatish Balay   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
841a7e14dcfSSatish Balay   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
842a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
843a7e14dcfSSatish Balay   PetscFunctionReturn(0);
844a7e14dcfSSatish Balay }
845a7e14dcfSSatish Balay 
846a7e14dcfSSatish Balay /*------------------------------------------------------------*/
8474416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptionItems *PetscOptionsObject,Tao tao)
848a7e14dcfSSatish Balay {
849a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
850a7e14dcfSSatish Balay   PetscErrorCode ierr;
851a7e14dcfSSatish Balay 
852a7e14dcfSSatish Balay   PetscFunctionBegin;
8531a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
854a7e14dcfSSatish 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);
855a7e14dcfSSatish 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);
856a7e14dcfSSatish 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);
857a7e14dcfSSatish 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);
85894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr);
85994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr);
86094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr);
86194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr);
86294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr);
86394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr);
86494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr);
86594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr);
86694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr);
86794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr);
86894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr);
86994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr);
87094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr);
87194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr);
87294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr);
87394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr);
87494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr);
87594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr);
87694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr);
87794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr);
87894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr);
87994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr);
88094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr);
88194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr);
88294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr);
88394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr);
88494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr);
88594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr);
88694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr);
88794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr);
88894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr);
88994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr);
89094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr);
89194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr);
89294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr);
89394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr);
89494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr);
89594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr);
89694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr);
89794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr);
89894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr);
89994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr);
90094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr);
90194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr);
90294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr);
903a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
904a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
905a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
906a7e14dcfSSatish Balay   PetscFunctionReturn(0);
907a7e14dcfSSatish Balay }
908a7e14dcfSSatish Balay 
909a7e14dcfSSatish Balay 
910a7e14dcfSSatish Balay /*------------------------------------------------------------*/
911441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
912a7e14dcfSSatish Balay {
913a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
914a7e14dcfSSatish Balay   PetscInt       nrejects;
915a7e14dcfSSatish Balay   PetscBool      isascii;
916a7e14dcfSSatish Balay   PetscErrorCode ierr;
917a7e14dcfSSatish Balay 
918a7e14dcfSSatish Balay   PetscFunctionBegin;
919a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
920a7e14dcfSSatish Balay   if (isascii) {
921a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
922a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
923*cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(nlsP->M,&nrejects);CHKERRQ(ierr);
924a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
925a7e14dcfSSatish Balay     }
926a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
927a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
928a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
929a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
930a7e14dcfSSatish Balay 
931a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
932a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
933a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
934a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
935a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
936a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
937a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
938a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
939a7e14dcfSSatish Balay   }
940a7e14dcfSSatish Balay   PetscFunctionReturn(0);
941a7e14dcfSSatish Balay }
942a7e14dcfSSatish Balay 
943a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
9444aa34175SJason Sarich /*MC
9454aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
9464aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
9474aa34175SJason Sarich   system of equations to obtain the step diretion dk:
9484aa34175SJason Sarich               Hk dk = -gk
9494aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
9504aa34175SJason Sarich   solve
9514aa34175SJason Sarich               min_t f(xk + t d_k)
9524aa34175SJason Sarich 
9534aa34175SJason Sarich     Options Database Keys:
9541daac38eSTodd Munson + -tao_nls_pc_type - "none","ahess","bfgs","petsc"
9554aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs"
9564aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation"
9574aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
9584aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
9594aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
9604aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
9614aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
9624aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
9634aa34175SJason Sarich . -tao_nls_pgfac - growth factor
9644aa34175SJason Sarich . -tao_nls_psfac - shrink factor
9654aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
9664aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
9674aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
9684aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
9694aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
9704aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
9714aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
9724aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
9734aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
9744aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
9754aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
9764aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
9774aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
9784aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
9794aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
9804aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
9814aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
9824aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
9834aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
9844aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
9854aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
9864aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
9874aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
9884aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
9891522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
9901522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
9911522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
9921522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
9931522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
9941522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
9951522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
9961eb8069cSJason Sarich 
9971eb8069cSJason Sarich   Level: beginner
9981522df2eSJason Sarich M*/
9994aa34175SJason Sarich 
1000728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
1001a7e14dcfSSatish Balay {
1002a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
10038caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
1004a7e14dcfSSatish Balay   PetscErrorCode ierr;
1005a7e14dcfSSatish Balay 
1006a7e14dcfSSatish Balay   PetscFunctionBegin;
10073c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1008a7e14dcfSSatish Balay 
1009a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
1010a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
1011a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
1012a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1013a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
1014a7e14dcfSSatish Balay 
10156552cf8aSJason Sarich   /* Override default settings (unless already changed) */
10166552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
10176552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
10186552cf8aSJason Sarich 
1019a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
1020a7e14dcfSSatish Balay 
1021a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
1022a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
1023a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
1024a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
1025a7e14dcfSSatish Balay 
1026a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
1027a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
1028a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
1029a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
1030a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
1031a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
1032a7e14dcfSSatish Balay 
1033a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
1034a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
1035a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
1036a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
1037a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
1038a7e14dcfSSatish Balay 
1039a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
1040a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
1041a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
1042a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
1043a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
1044a7e14dcfSSatish Balay 
1045a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
1046a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
1047a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
1048a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
1049a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
1050a7e14dcfSSatish Balay 
1051a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
1052a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
1053a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
1054a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
1055a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
1056a7e14dcfSSatish Balay 
1057a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
1058a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
1059a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
1060a7e14dcfSSatish Balay 
1061a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
1062a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
1063a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
1064a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
1065a7e14dcfSSatish Balay 
1066a7e14dcfSSatish Balay   nlsP->theta = 0.05;
1067a7e14dcfSSatish Balay 
1068a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
1069a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
1070a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
1071a7e14dcfSSatish Balay 
1072a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
1073a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
1074a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
1075a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
1076a7e14dcfSSatish Balay 
1077a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
1078a7e14dcfSSatish Balay 
1079a7e14dcfSSatish Balay   /*  Remaining parameters */
1080a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
1081a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
1082a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
1083a7e14dcfSSatish Balay 
1084a7e14dcfSSatish Balay   nlsP->pc_type         = NLS_PC_BFGS;
1085a7e14dcfSSatish Balay   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1086a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1087a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
1088a7e14dcfSSatish Balay 
1089a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
109063b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1091a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1092441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
10935d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1094a7e14dcfSSatish Balay 
1095a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
1096a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
109763b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
10985d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
10991daac38eSTodd Munson   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1100a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1101a7e14dcfSSatish Balay }
1102