xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision aaa7dc30da3270cff6cb10b1db605b2ca746f216)
1*aaa7dc30SBarry Smith #include <taolinesearch.h>
2*aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
3*aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/nls/nls.h>
4a7e14dcfSSatish Balay 
5*aaa7dc30SBarry Smith #include <petscksp.h>
6*aaa7dc30SBarry Smith #include <petscpc.h>
7*aaa7dc30SBarry Smith #include <petsc-private/kspimpl.h>
8*aaa7dc30SBarry Smith #include <petsc-private/pcimpl.h>
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay #define NLS_KSP_CG      0
11a7e14dcfSSatish Balay #define NLS_KSP_NASH    1
12a7e14dcfSSatish Balay #define NLS_KSP_STCG    2
13a7e14dcfSSatish Balay #define NLS_KSP_GLTR    3
14a7e14dcfSSatish Balay #define NLS_KSP_PETSC   4
15a7e14dcfSSatish Balay #define NLS_KSP_TYPES   5
16a7e14dcfSSatish Balay 
17a7e14dcfSSatish Balay #define NLS_PC_NONE     0
18a7e14dcfSSatish Balay #define NLS_PC_AHESS    1
19a7e14dcfSSatish Balay #define NLS_PC_BFGS     2
20a7e14dcfSSatish Balay #define NLS_PC_PETSC    3
21a7e14dcfSSatish Balay #define NLS_PC_TYPES    4
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS        0
24a7e14dcfSSatish Balay #define BFGS_SCALE_PHESS        1
25a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS         2
26a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES        3
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT         0
29a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION        1
30a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION    2
31a7e14dcfSSatish Balay #define NLS_INIT_TYPES            3
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay #define NLS_UPDATE_STEP           0
34a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION      1
35a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION  2
36a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES          3
37a7e14dcfSSatish Balay 
3887f595a5SBarry Smith static const char *NLS_KSP[64] = {"cg", "nash", "stcg", "gltr", "petsc"};
39a7e14dcfSSatish Balay 
4087f595a5SBarry Smith static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"};
41a7e14dcfSSatish Balay 
4287f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
43a7e14dcfSSatish Balay 
4487f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"};
45a7e14dcfSSatish Balay 
4687f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"};
47a7e14dcfSSatish Balay 
48a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x);
49a7e14dcfSSatish Balay /* Routine for BFGS preconditioner
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay  Implements Newton's Method with a line search approach for solving
53a7e14dcfSSatish Balay  unconstrained minimization problems.  A More'-Thuente line search
54a7e14dcfSSatish Balay  is used to guarantee that the bfgs preconditioner remains positive
55a7e14dcfSSatish Balay  definite.
56a7e14dcfSSatish Balay 
57a7e14dcfSSatish Balay  The method can shift the Hessian matrix.  The shifting procedure is
58a7e14dcfSSatish Balay  adapted from the PATH algorithm for solving complementarity
59a7e14dcfSSatish Balay  problems.
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay  The linear system solve should be done with a conjugate gradient
62a7e14dcfSSatish Balay  method, although any method can be used. */
63a7e14dcfSSatish Balay 
64a7e14dcfSSatish Balay #define NLS_NEWTON              0
65a7e14dcfSSatish Balay #define NLS_BFGS                1
66a7e14dcfSSatish Balay #define NLS_SCALED_GRADIENT     2
67a7e14dcfSSatish Balay #define NLS_GRADIENT            3
68a7e14dcfSSatish Balay 
69a7e14dcfSSatish Balay #undef __FUNCT__
70a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NLS"
71a7e14dcfSSatish Balay static PetscErrorCode TaoSolve_NLS(TaoSolver tao)
72a7e14dcfSSatish Balay {
73a7e14dcfSSatish Balay   PetscErrorCode                 ierr;
74a7e14dcfSSatish Balay   TAO_NLS                        *nlsP = (TAO_NLS *)tao->data;
75a7e14dcfSSatish Balay   PC                             pc;
76a7e14dcfSSatish Balay   KSPConvergedReason             ksp_reason;
77a7e14dcfSSatish Balay   TaoLineSearchTerminationReason ls_reason;
78a7e14dcfSSatish Balay   TaoSolverTerminationReason     reason;
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay   PetscReal                      fmin, ftrial, f_full, prered, actred, kappa, sigma;
81a7e14dcfSSatish Balay   PetscReal                      tau, tau_1, tau_2, tau_max, tau_min, max_radius;
82a7e14dcfSSatish Balay   PetscReal                      f, fold, gdx, gnorm, pert;
83a7e14dcfSSatish Balay   PetscReal                      step = 1.0;
84a7e14dcfSSatish Balay   PetscReal                      delta;
85a7e14dcfSSatish Balay   PetscReal                      norm_d = 0.0, e_min;
86a7e14dcfSSatish Balay 
87a7e14dcfSSatish Balay   MatStructure                   matflag;
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay   PetscInt                       stepType;
90a7e14dcfSSatish Balay   PetscInt                       iter = 0;
91a7e14dcfSSatish Balay   PetscInt                       bfgsUpdates = 0;
92a7e14dcfSSatish Balay   PetscInt                       n,N,kspits;
93a7e14dcfSSatish Balay   PetscInt                       needH;
94a7e14dcfSSatish Balay 
95a7e14dcfSSatish Balay   PetscInt                       i_max = 5;
96a7e14dcfSSatish Balay   PetscInt                       j_max = 1;
97a7e14dcfSSatish Balay   PetscInt                       i, j;
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay   PetscFunctionBegin;
100a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
101a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
102a7e14dcfSSatish Balay   }
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay   /* Initialized variables */
105a7e14dcfSSatish Balay   pert = nlsP->sval;
106a7e14dcfSSatish Balay 
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   /* Modify the linear solver to a trust region method if desired */
116a7e14dcfSSatish Balay   switch(nlsP->ksp_type) {
117a7e14dcfSSatish Balay   case NLS_KSP_CG:
118a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr);
119a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
120a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
121a7e14dcfSSatish Balay     }
122a7e14dcfSSatish Balay     break;
123a7e14dcfSSatish Balay 
124a7e14dcfSSatish Balay   case NLS_KSP_NASH:
125a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
126a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
127a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
128a7e14dcfSSatish Balay     }
129a7e14dcfSSatish Balay     break;
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay   case NLS_KSP_STCG:
132a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
133a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
134a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
135a7e14dcfSSatish Balay     }
136a7e14dcfSSatish Balay     break;
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay   case NLS_KSP_GLTR:
139a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
140a7e14dcfSSatish Balay     if (tao->ksp->ops->setfromoptions) {
141a7e14dcfSSatish Balay       (*tao->ksp->ops->setfromoptions)(tao->ksp);
142a7e14dcfSSatish Balay     }
143a7e14dcfSSatish Balay     break;
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   default:
146a7e14dcfSSatish Balay     /* Use the method set by the ksp_type */
147a7e14dcfSSatish Balay     break;
148a7e14dcfSSatish Balay   }
149a7e14dcfSSatish Balay 
150a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
151a7e14dcfSSatish Balay    Will be reset during the first iteration */
152a7e14dcfSSatish Balay   if (NLS_KSP_NASH == nlsP->ksp_type) {
153a7e14dcfSSatish Balay     ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   } else if (NLS_KSP_STCG == nlsP->ksp_type) {
155a7e14dcfSSatish Balay     ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
156a7e14dcfSSatish Balay   } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
157a7e14dcfSSatish Balay     ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
158a7e14dcfSSatish Balay   }
159a7e14dcfSSatish Balay 
16087f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
161a7e14dcfSSatish Balay     tao->trust = tao->trust0;
162a7e14dcfSSatish Balay 
16387f595a5SBarry Smith     if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative");
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay     /* Modify the radius if it is too large or small */
166a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
167a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
168a7e14dcfSSatish Balay   }
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay   /* Get vectors we will need */
171a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
172a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
173a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
174a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
175a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr);
176a7e14dcfSSatish Balay   }
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay   /* Check convergence criteria */
179a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
180a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
18187f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
182a7e14dcfSSatish Balay   needH = 1;
183a7e14dcfSSatish Balay 
184a7e14dcfSSatish Balay   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
18587f595a5SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay   /* create vectors for the limited memory preconditioner */
18887f595a5SBarry Smith   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
189a7e14dcfSSatish Balay     if (!nlsP->Diag) {
190a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
191a7e14dcfSSatish Balay     }
192a7e14dcfSSatish Balay   }
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay   /* Modify the preconditioner to use the bfgs approximation */
195a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
196a7e14dcfSSatish Balay   switch(nlsP->pc_type) {
197a7e14dcfSSatish Balay   case NLS_PC_NONE:
198a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
199a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
200a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
201a7e14dcfSSatish Balay     }
202a7e14dcfSSatish Balay     break;
203a7e14dcfSSatish Balay 
204a7e14dcfSSatish Balay   case NLS_PC_AHESS:
205a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
206a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
207a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
208a7e14dcfSSatish Balay     }
209a7e14dcfSSatish Balay     ierr = PCJacobiSetUseAbs(pc);CHKERRQ(ierr);
210a7e14dcfSSatish Balay     break;
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay   case NLS_PC_BFGS:
213a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
214a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
215a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
216a7e14dcfSSatish Balay     }
217a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
218a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
219a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
220a7e14dcfSSatish Balay     break;
221a7e14dcfSSatish Balay 
222a7e14dcfSSatish Balay   default:
223a7e14dcfSSatish Balay     /* Use the pc method set by pc_type */
224a7e14dcfSSatish Balay     break;
225a7e14dcfSSatish Balay   }
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
228a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
22987f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
230a7e14dcfSSatish Balay     switch(nlsP->init_type) {
231a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
232a7e14dcfSSatish Balay       /* Use the initial radius specified */
233a7e14dcfSSatish Balay       break;
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
236a7e14dcfSSatish Balay       /* Use the initial radius specified */
237a7e14dcfSSatish Balay       max_radius = 0.0;
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
240a7e14dcfSSatish Balay         fmin = f;
241a7e14dcfSSatish Balay         sigma = 0.0;
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay         if (needH) {
244a7e14dcfSSatish Balay           ierr = TaoComputeHessian(tao, tao->solution, &tao->hessian, &tao->hessian_pre, &matflag);CHKERRQ(ierr);
245a7e14dcfSSatish Balay           needH = 0;
246a7e14dcfSSatish Balay         }
247a7e14dcfSSatish Balay 
248a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
249a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
250a7e14dcfSSatish Balay           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
251a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
252a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
253a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
25487f595a5SBarry Smith           } else {
255a7e14dcfSSatish Balay             if (ftrial < fmin) {
256a7e14dcfSSatish Balay               fmin = ftrial;
257a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
258a7e14dcfSSatish Balay             }
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
261a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
262a7e14dcfSSatish Balay 
263a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
264a7e14dcfSSatish Balay             actred = f - ftrial;
26587f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
266a7e14dcfSSatish Balay               kappa = 1.0;
26787f595a5SBarry Smith             } else {
268a7e14dcfSSatish Balay               kappa = actred / prered;
269a7e14dcfSSatish Balay             }
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
272a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
273a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
274a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
275a7e14dcfSSatish Balay 
276a7e14dcfSSatish Balay             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
277a7e14dcfSSatish Balay               /* Great agreement */
278a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
279a7e14dcfSSatish Balay 
280a7e14dcfSSatish Balay               if (tau_max < 1.0) {
281a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
28287f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
283a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
28487f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
285a7e14dcfSSatish Balay                 tau = tau_1;
28687f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
287a7e14dcfSSatish Balay                 tau = tau_2;
28887f595a5SBarry Smith               } else {
289a7e14dcfSSatish Balay                 tau = tau_max;
290a7e14dcfSSatish Balay               }
29187f595a5SBarry Smith             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
292a7e14dcfSSatish Balay               /* Good agreement */
293a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
296a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
29787f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
298a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
29987f595a5SBarry Smith               } else {
300a7e14dcfSSatish Balay                 tau = tau_max;
301a7e14dcfSSatish Balay               }
30287f595a5SBarry Smith             } else {
303a7e14dcfSSatish Balay               /* Not good agreement */
304a7e14dcfSSatish Balay               if (tau_min > 1.0) {
305a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
30687f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
307a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
30887f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
309a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
31087f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
311a7e14dcfSSatish Balay                 tau = tau_1;
31287f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
313a7e14dcfSSatish Balay                 tau = tau_2;
31487f595a5SBarry Smith               } else {
315a7e14dcfSSatish Balay                 tau = tau_max;
316a7e14dcfSSatish Balay               }
317a7e14dcfSSatish Balay             }
318a7e14dcfSSatish Balay           }
319a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
320a7e14dcfSSatish Balay         }
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay         if (fmin < f) {
323a7e14dcfSSatish Balay           f = fmin;
324a7e14dcfSSatish Balay           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
325a7e14dcfSSatish Balay           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
326a7e14dcfSSatish Balay 
327a7e14dcfSSatish Balay           ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
32887f595a5SBarry Smith           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
329a7e14dcfSSatish Balay           needH = 1;
330a7e14dcfSSatish Balay 
331a7e14dcfSSatish Balay           ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
33287f595a5SBarry Smith           if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
333a7e14dcfSSatish Balay         }
334a7e14dcfSSatish Balay       }
335a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
336a7e14dcfSSatish Balay 
337a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
338a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
339a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
340a7e14dcfSSatish Balay       break;
341a7e14dcfSSatish Balay 
342a7e14dcfSSatish Balay     default:
343a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
344a7e14dcfSSatish Balay       tao->trust = 0.0;
345a7e14dcfSSatish Balay       break;
346a7e14dcfSSatish Balay     }
347a7e14dcfSSatish Balay   }
348a7e14dcfSSatish Balay 
349a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
350a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
351a7e14dcfSSatish Balay      since the function value may have decreased */
352a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type) {
353a7e14dcfSSatish Balay     if (f != 0.0) {
354a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
35587f595a5SBarry Smith     } else {
356a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
357a7e14dcfSSatish Balay     }
358a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
359a7e14dcfSSatish Balay   }
360a7e14dcfSSatish Balay 
361a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
362a7e14dcfSSatish Balay   nlsP->newt = 0;
363a7e14dcfSSatish Balay   nlsP->bfgs = 0;
364a7e14dcfSSatish Balay   nlsP->sgrad = 0;
365a7e14dcfSSatish Balay   nlsP->grad = 0;
366a7e14dcfSSatish Balay 
367a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
368a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
369a7e14dcfSSatish Balay     ++iter;
370a7e14dcfSSatish Balay 
371a7e14dcfSSatish Balay     /* Compute the Hessian */
372a7e14dcfSSatish Balay     if (needH) {
373a7e14dcfSSatish Balay       ierr = TaoComputeHessian(tao, tao->solution, &tao->hessian, &tao->hessian_pre, &matflag);CHKERRQ(ierr);
374a7e14dcfSSatish Balay       needH = 0;
375a7e14dcfSSatish Balay     }
376a7e14dcfSSatish Balay 
37787f595a5SBarry Smith     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
378a7e14dcfSSatish Balay       /* Obtain diagonal for the bfgs preconditioner  */
379a7e14dcfSSatish Balay       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
380a7e14dcfSSatish Balay       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
381a7e14dcfSSatish Balay       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
382a7e14dcfSSatish Balay       ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
383a7e14dcfSSatish Balay     }
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
386a7e14dcfSSatish Balay     if (pert > 0) {
387a7e14dcfSSatish Balay       ierr = MatShift(tao->hessian, pert);
388a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
389a7e14dcfSSatish Balay         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
390a7e14dcfSSatish Balay       }
391a7e14dcfSSatish Balay     }
392a7e14dcfSSatish Balay 
393a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type) {
394a7e14dcfSSatish Balay       if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) {
395a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner  */
396a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
397a7e14dcfSSatish Balay         ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
398a7e14dcfSSatish Balay         ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
399a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
400a7e14dcfSSatish Balay       }
401a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
402a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
403a7e14dcfSSatish Balay       ++bfgsUpdates;
404a7e14dcfSSatish Balay     }
405a7e14dcfSSatish Balay 
406a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
407a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre,matflag);CHKERRQ(ierr);
40887f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type ||  NLS_KSP_GLTR == nlsP->ksp_type) {
409a7e14dcfSSatish Balay 
410a7e14dcfSSatish Balay       if (NLS_KSP_NASH == nlsP->ksp_type) {
411a7e14dcfSSatish Balay         ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
412a7e14dcfSSatish Balay       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
413a7e14dcfSSatish Balay          ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
414a7e14dcfSSatish Balay       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
415a7e14dcfSSatish Balay         ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
416a7e14dcfSSatish Balay       }
417a7e14dcfSSatish Balay 
418a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
419a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
420a7e14dcfSSatish Balay       tao->ksp_its+=kspits;
421a7e14dcfSSatish Balay 
422a7e14dcfSSatish Balay       if (NLS_KSP_NASH == nlsP->ksp_type) {
423a7e14dcfSSatish Balay         ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
424a7e14dcfSSatish Balay       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
425a7e14dcfSSatish Balay          ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
426a7e14dcfSSatish Balay       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
427a7e14dcfSSatish Balay         ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
428a7e14dcfSSatish Balay       }
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
431a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
432a7e14dcfSSatish Balay         if (norm_d > 0.0) {
433a7e14dcfSSatish Balay           tao->trust = norm_d;
434a7e14dcfSSatish Balay 
435a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
436a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
437a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
43887f595a5SBarry Smith         } else {
439a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
440a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
441a7e14dcfSSatish Balay           tao->trust = tao->trust0;
442a7e14dcfSSatish Balay 
443a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
444a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
445a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
446a7e14dcfSSatish Balay 
447a7e14dcfSSatish Balay           if (NLS_KSP_NASH == nlsP->ksp_type) {
448a7e14dcfSSatish Balay             ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
449a7e14dcfSSatish Balay           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
450a7e14dcfSSatish Balay             ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
451a7e14dcfSSatish Balay           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
452a7e14dcfSSatish Balay             ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
453a7e14dcfSSatish Balay           }
454a7e14dcfSSatish Balay 
455a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
456a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
457a7e14dcfSSatish Balay           tao->ksp_its+=kspits;
458a7e14dcfSSatish Balay           if (NLS_KSP_NASH == nlsP->ksp_type) {
459a7e14dcfSSatish Balay             ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
460a7e14dcfSSatish Balay           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
461a7e14dcfSSatish Balay             ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
462a7e14dcfSSatish Balay           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
463a7e14dcfSSatish Balay             ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
464a7e14dcfSSatish Balay           }
465a7e14dcfSSatish Balay 
46687f595a5SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
467a7e14dcfSSatish Balay         }
468a7e14dcfSSatish Balay       }
46987f595a5SBarry Smith     } else {
470a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
471a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
472a7e14dcfSSatish Balay       tao->ksp_its += kspits;
473a7e14dcfSSatish Balay     }
474a7e14dcfSSatish Balay     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
475a7e14dcfSSatish Balay     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
47687f595a5SBarry Smith     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
477a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
478a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
479a7e14dcfSSatish Balay 
480a7e14dcfSSatish Balay       if (f != 0.0) {
481a7e14dcfSSatish Balay         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
48287f595a5SBarry Smith       } else {
483a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
484a7e14dcfSSatish Balay       }
485a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
486a7e14dcfSSatish Balay       ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
487a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
488a7e14dcfSSatish Balay       bfgsUpdates = 1;
489a7e14dcfSSatish Balay     }
490a7e14dcfSSatish Balay 
491a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
492a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
49387f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
494a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
49587f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
496a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
49787f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
498a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
49987f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
500a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
50187f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
502a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
50387f595a5SBarry Smith     } else {
504a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
505a7e14dcfSSatish Balay     }
506a7e14dcfSSatish Balay 
507a7e14dcfSSatish Balay     /* Check for success (descent direction) */
508a7e14dcfSSatish Balay     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
509a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
510a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
511a7e14dcfSSatish Balay          Update the perturbation for next time */
512a7e14dcfSSatish Balay       if (pert <= 0.0) {
513a7e14dcfSSatish Balay         /* Initialize the perturbation */
514a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
515a7e14dcfSSatish Balay         if (NLS_KSP_GLTR == nlsP->ksp_type) {
516a7e14dcfSSatish Balay           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
517a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
518a7e14dcfSSatish Balay         }
51987f595a5SBarry Smith       } else {
520a7e14dcfSSatish Balay         /* Increase the perturbation */
521a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
522a7e14dcfSSatish Balay       }
523a7e14dcfSSatish Balay 
524a7e14dcfSSatish Balay       if (NLS_PC_BFGS != nlsP->pc_type) {
525a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
526a7e14dcfSSatish Balay            Must use gradient direction in this case */
527a7e14dcfSSatish Balay         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
528a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
529a7e14dcfSSatish Balay         ++nlsP->grad;
530a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
53187f595a5SBarry Smith       } else {
532a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
533a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
534a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
535a7e14dcfSSatish Balay 
536a7e14dcfSSatish Balay         /* Check for success (descent direction) */
537a7e14dcfSSatish Balay         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
538a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
539a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
540a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
541a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
542a7e14dcfSSatish Balay              which is guaranteed to be descent */
543a7e14dcfSSatish Balay 
544a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
545a7e14dcfSSatish Balay 
546a7e14dcfSSatish Balay           if (f != 0.0) {
547a7e14dcfSSatish Balay             delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
54887f595a5SBarry Smith           } else {
549a7e14dcfSSatish Balay             delta = 2.0 / (gnorm*gnorm);
550a7e14dcfSSatish Balay           }
551a7e14dcfSSatish Balay           ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
552a7e14dcfSSatish Balay           ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
553a7e14dcfSSatish Balay           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
554a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
555a7e14dcfSSatish Balay           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
556a7e14dcfSSatish Balay 
557a7e14dcfSSatish Balay           bfgsUpdates = 1;
558a7e14dcfSSatish Balay           ++nlsP->sgrad;
559a7e14dcfSSatish Balay           stepType = NLS_SCALED_GRADIENT;
56087f595a5SBarry Smith         } else {
561a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
562a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
563a7e14dcfSSatish Balay             ++nlsP->sgrad;
564a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
56587f595a5SBarry Smith           } else {
566a7e14dcfSSatish Balay             ++nlsP->bfgs;
567a7e14dcfSSatish Balay             stepType = NLS_BFGS;
568a7e14dcfSSatish Balay           }
569a7e14dcfSSatish Balay         }
570a7e14dcfSSatish Balay       }
57187f595a5SBarry Smith     } else {
572a7e14dcfSSatish Balay       /* Computed Newton step is descent */
573a7e14dcfSSatish Balay       switch (ksp_reason) {
574a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
575a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
576a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
577a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
578a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
579a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
580a7e14dcfSSatish Balay         if (pert <= 0.0) {
581a7e14dcfSSatish Balay           /* Initialize the perturbation */
582a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
583a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
584a7e14dcfSSatish Balay             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
585a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
586a7e14dcfSSatish Balay           }
58787f595a5SBarry Smith         } else {
588a7e14dcfSSatish Balay           /* Increase the perturbation */
589a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
590a7e14dcfSSatish Balay         }
591a7e14dcfSSatish Balay         break;
592a7e14dcfSSatish Balay 
593a7e14dcfSSatish Balay       default:
594a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
595a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
596a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
597a7e14dcfSSatish Balay           pert = 0.0;
598a7e14dcfSSatish Balay         }
599a7e14dcfSSatish Balay         break;
600a7e14dcfSSatish Balay       }
601a7e14dcfSSatish Balay 
602a7e14dcfSSatish Balay       ++nlsP->newt;
603a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
604a7e14dcfSSatish Balay     }
605a7e14dcfSSatish Balay 
606a7e14dcfSSatish Balay     /* Perform the linesearch */
607a7e14dcfSSatish Balay     fold = f;
608a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
609a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
610a7e14dcfSSatish Balay 
611a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
612a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
613a7e14dcfSSatish Balay 
61487f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
615a7e14dcfSSatish Balay       f = fold;
616a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
617a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
618a7e14dcfSSatish Balay 
619a7e14dcfSSatish Balay       switch(stepType) {
620a7e14dcfSSatish Balay       case NLS_NEWTON:
621a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
622a7e14dcfSSatish Balay            Update the perturbation for next time */
623a7e14dcfSSatish Balay         if (pert <= 0.0) {
624a7e14dcfSSatish Balay           /* Initialize the perturbation */
625a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
626a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
627a7e14dcfSSatish Balay             ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
628a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
629a7e14dcfSSatish Balay           }
63087f595a5SBarry Smith         } else {
631a7e14dcfSSatish Balay           /* Increase the perturbation */
632a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
633a7e14dcfSSatish Balay         }
634a7e14dcfSSatish Balay 
635a7e14dcfSSatish Balay         if (NLS_PC_BFGS != nlsP->pc_type) {
636a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
637a7e14dcfSSatish Balay              Must use gradient direction in this case */
638a7e14dcfSSatish Balay           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
639a7e14dcfSSatish Balay           ++nlsP->grad;
640a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
64187f595a5SBarry Smith         } else {
642a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
643a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
644a7e14dcfSSatish Balay           /* Check for success (descent direction) */
645a7e14dcfSSatish Balay           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
646a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
647a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
648a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
649a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
650a7e14dcfSSatish Balay 
651a7e14dcfSSatish Balay             if (f != 0.0) {
652a7e14dcfSSatish Balay               delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
65387f595a5SBarry Smith             } else {
654a7e14dcfSSatish Balay               delta = 2.0 / (gnorm*gnorm);
655a7e14dcfSSatish Balay             }
656a7e14dcfSSatish Balay             ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
657a7e14dcfSSatish Balay             ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
658a7e14dcfSSatish Balay             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
659a7e14dcfSSatish Balay             ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
660a7e14dcfSSatish Balay 
661a7e14dcfSSatish Balay             bfgsUpdates = 1;
662a7e14dcfSSatish Balay             ++nlsP->sgrad;
663a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
66487f595a5SBarry Smith           } else {
665a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
666a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
667a7e14dcfSSatish Balay               ++nlsP->sgrad;
668a7e14dcfSSatish Balay               stepType = NLS_SCALED_GRADIENT;
66987f595a5SBarry Smith             } else {
670a7e14dcfSSatish Balay               ++nlsP->bfgs;
671a7e14dcfSSatish Balay               stepType = NLS_BFGS;
672a7e14dcfSSatish Balay             }
673a7e14dcfSSatish Balay           }
674a7e14dcfSSatish Balay         }
675a7e14dcfSSatish Balay         break;
676a7e14dcfSSatish Balay 
677a7e14dcfSSatish Balay       case NLS_BFGS:
678a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
679a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
680a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
681a7e14dcfSSatish Balay 
682a7e14dcfSSatish Balay         if (f != 0.0) {
683a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
68487f595a5SBarry Smith         } else {
685a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
686a7e14dcfSSatish Balay         }
687a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
688a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
689a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
690a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
691a7e14dcfSSatish Balay 
692a7e14dcfSSatish Balay         bfgsUpdates = 1;
693a7e14dcfSSatish Balay         ++nlsP->sgrad;
694a7e14dcfSSatish Balay         stepType = NLS_SCALED_GRADIENT;
695a7e14dcfSSatish Balay         break;
696a7e14dcfSSatish Balay 
697a7e14dcfSSatish Balay       case NLS_SCALED_GRADIENT:
698a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
699a7e14dcfSSatish Balay            The scaled gradient step did not produce a new iterate;
700a7e14dcfSSatish Balay            attemp to use the gradient direction.
701a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
702a7e14dcfSSatish Balay 
703a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr);
704a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr);
705a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
706a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
707a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
708a7e14dcfSSatish Balay 
709a7e14dcfSSatish Balay         bfgsUpdates = 1;
710a7e14dcfSSatish Balay         ++nlsP->grad;
711a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
712a7e14dcfSSatish Balay         break;
713a7e14dcfSSatish Balay       }
714a7e14dcfSSatish Balay       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
715a7e14dcfSSatish Balay 
716a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
717a7e14dcfSSatish Balay       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
718a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
719a7e14dcfSSatish Balay     }
720a7e14dcfSSatish Balay 
72187f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
722a7e14dcfSSatish Balay       /* Failed to find an improving point */
723a7e14dcfSSatish Balay       f = fold;
724a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
725a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
726a7e14dcfSSatish Balay       step = 0.0;
727a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
728a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
729a7e14dcfSSatish Balay       break;
730a7e14dcfSSatish Balay     }
731a7e14dcfSSatish Balay 
732a7e14dcfSSatish Balay     /* Update trust region radius */
73387f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
734a7e14dcfSSatish Balay       switch(nlsP->update_type) {
735a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
736a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
737a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
738a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
739a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
74087f595a5SBarry Smith           } else if (step < nlsP->nu2) {
741a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
742a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
74387f595a5SBarry Smith           } else if (step < nlsP->nu3) {
744a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
745a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
746a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
74787f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
748a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
749a7e14dcfSSatish Balay             }
75087f595a5SBarry Smith           } else if (step < nlsP->nu4) {
751a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
752a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
75387f595a5SBarry Smith           } else {
754a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
755a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
756a7e14dcfSSatish Balay           }
75787f595a5SBarry Smith         } else {
758a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
759a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
760a7e14dcfSSatish Balay         }
761a7e14dcfSSatish Balay         break;
762a7e14dcfSSatish Balay 
763a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
764a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
765a7e14dcfSSatish Balay           /*  Get predicted reduction */
766a7e14dcfSSatish Balay 
767a7e14dcfSSatish Balay           if (NLS_KSP_STCG == nlsP->ksp_type) {
768a7e14dcfSSatish Balay             ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);
769a7e14dcfSSatish Balay           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
770a7e14dcfSSatish Balay             ierr = KSPNASHGetObjFcn(tao->ksp,&prered);
771a7e14dcfSSatish Balay           } else {
772a7e14dcfSSatish Balay             ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
773a7e14dcfSSatish Balay           }
774a7e14dcfSSatish Balay 
775a7e14dcfSSatish Balay           if (prered >= 0.0) {
776a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
777a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
778a7e14dcfSSatish Balay             /*  be rejected! */
779a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
78087f595a5SBarry Smith           } else {
781a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
782a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
78387f595a5SBarry Smith             } else {
784a7e14dcfSSatish Balay               /*  Compute and actual reduction */
785a7e14dcfSSatish Balay               actred = fold - f_full;
786a7e14dcfSSatish Balay               prered = -prered;
787a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
788a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
789a7e14dcfSSatish Balay                 kappa = 1.0;
79087f595a5SBarry Smith               } else {
791a7e14dcfSSatish Balay                 kappa = actred / prered;
792a7e14dcfSSatish Balay               }
793a7e14dcfSSatish Balay 
794a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
795a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
796a7e14dcfSSatish Balay                 /*  Very bad step */
797a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
79887f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
799a7e14dcfSSatish Balay                 /*  Marginal bad step */
800a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
80187f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
802a7e14dcfSSatish Balay                 /*  Reasonable step */
803a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
804a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
80587f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
806a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
807a7e14dcfSSatish Balay                 }
80887f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
809a7e14dcfSSatish Balay                 /*  Good step */
810a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
81187f595a5SBarry Smith               } else {
812a7e14dcfSSatish Balay                 /*  Very good step */
813a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
814a7e14dcfSSatish Balay               }
815a7e14dcfSSatish Balay             }
816a7e14dcfSSatish Balay           }
81787f595a5SBarry Smith         } else {
818a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
819a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
820a7e14dcfSSatish Balay         }
821a7e14dcfSSatish Balay         break;
822a7e14dcfSSatish Balay 
823a7e14dcfSSatish Balay       default:
824a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
825a7e14dcfSSatish Balay 
826a7e14dcfSSatish Balay           if (NLS_KSP_STCG == nlsP->ksp_type) {
827a7e14dcfSSatish Balay               ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);
828a7e14dcfSSatish Balay           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
829a7e14dcfSSatish Balay               ierr = KSPNASHGetObjFcn(tao->ksp,&prered);
830a7e14dcfSSatish Balay           } else {
831a7e14dcfSSatish Balay               ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
832a7e14dcfSSatish Balay           }
833a7e14dcfSSatish Balay           if (prered >= 0.0) {
834a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
835a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
836a7e14dcfSSatish Balay             /*  be rejected! */
837a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
83887f595a5SBarry Smith           } else {
839a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
840a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
84187f595a5SBarry Smith             } else {
842a7e14dcfSSatish Balay               actred = fold - f_full;
843a7e14dcfSSatish Balay               prered = -prered;
84487f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
845a7e14dcfSSatish Balay                 kappa = 1.0;
84687f595a5SBarry Smith               } else {
847a7e14dcfSSatish Balay                 kappa = actred / prered;
848a7e14dcfSSatish Balay               }
849a7e14dcfSSatish Balay 
850a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
851a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
852a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
853a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
854a7e14dcfSSatish Balay 
855a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
856a7e14dcfSSatish Balay                 /*  Great agreement */
857a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
858a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
85987f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
860a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
86187f595a5SBarry Smith                 } else {
862a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
863a7e14dcfSSatish Balay                 }
86487f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
865a7e14dcfSSatish Balay                 /*  Good agreement */
866a7e14dcfSSatish Balay 
867a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
868a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
86987f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
870a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
87187f595a5SBarry Smith                 } else if (tau_max < 1.0) {
872a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
87387f595a5SBarry Smith                 } else {
874a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
875a7e14dcfSSatish Balay                 }
87687f595a5SBarry Smith               } else {
877a7e14dcfSSatish Balay                 /*  Not good agreement */
878a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
879a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
88087f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
881a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
88287f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
883a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
88487f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
885a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
88687f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
887a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
88887f595a5SBarry Smith                 } else {
889a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
890a7e14dcfSSatish Balay                 }
891a7e14dcfSSatish Balay               }
892a7e14dcfSSatish Balay             }
893a7e14dcfSSatish Balay           }
89487f595a5SBarry Smith         } else {
895a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
896a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
897a7e14dcfSSatish Balay         }
898a7e14dcfSSatish Balay         break;
899a7e14dcfSSatish Balay       }
900a7e14dcfSSatish Balay 
901a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
902a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
903a7e14dcfSSatish Balay     }
904a7e14dcfSSatish Balay 
905a7e14dcfSSatish Balay     /*  Check for termination */
906a7e14dcfSSatish Balay     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
90787f595a5SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
908a7e14dcfSSatish Balay     needH = 1;
909a7e14dcfSSatish Balay     ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
910a7e14dcfSSatish Balay   }
911a7e14dcfSSatish Balay   PetscFunctionReturn(0);
912a7e14dcfSSatish Balay }
913a7e14dcfSSatish Balay 
914a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
915a7e14dcfSSatish Balay #undef __FUNCT__
916a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NLS"
917a7e14dcfSSatish Balay static PetscErrorCode TaoSetUp_NLS(TaoSolver tao)
918a7e14dcfSSatish Balay {
919a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
920a7e14dcfSSatish Balay   PetscErrorCode ierr;
921a7e14dcfSSatish Balay 
922a7e14dcfSSatish Balay   PetscFunctionBegin;
923a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
924a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
925a7e14dcfSSatish Balay   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
926a7e14dcfSSatish Balay   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
927a7e14dcfSSatish Balay   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
928a7e14dcfSSatish Balay   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
929a7e14dcfSSatish Balay   nlsP->Diag = 0;
930a7e14dcfSSatish Balay   nlsP->M = 0;
931a7e14dcfSSatish Balay   PetscFunctionReturn(0);
932a7e14dcfSSatish Balay }
933a7e14dcfSSatish Balay 
934a7e14dcfSSatish Balay /*------------------------------------------------------------*/
935a7e14dcfSSatish Balay #undef __FUNCT__
936a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NLS"
937a7e14dcfSSatish Balay static PetscErrorCode TaoDestroy_NLS(TaoSolver tao)
938a7e14dcfSSatish Balay {
939a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
940a7e14dcfSSatish Balay   PetscErrorCode ierr;
941a7e14dcfSSatish Balay 
942a7e14dcfSSatish Balay   PetscFunctionBegin;
943a7e14dcfSSatish Balay   if (tao->setupcalled) {
944a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
945a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
946a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
947a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
948a7e14dcfSSatish Balay   }
949a7e14dcfSSatish Balay   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
950a7e14dcfSSatish Balay   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
951a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
952a7e14dcfSSatish Balay   PetscFunctionReturn(0);
953a7e14dcfSSatish Balay }
954a7e14dcfSSatish Balay 
955a7e14dcfSSatish Balay /*------------------------------------------------------------*/
956a7e14dcfSSatish Balay #undef __FUNCT__
957a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NLS"
958a7e14dcfSSatish Balay static PetscErrorCode TaoSetFromOptions_NLS(TaoSolver tao)
959a7e14dcfSSatish Balay {
960a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
961a7e14dcfSSatish Balay   PetscErrorCode ierr;
962a7e14dcfSSatish Balay 
963a7e14dcfSSatish Balay   PetscFunctionBegin;
964a7e14dcfSSatish Balay   ierr = PetscOptionsHead("Newton line search method for unconstrained optimization");CHKERRQ(ierr);
965a7e14dcfSSatish Balay   ierr = PetscOptionsEList("-tao_nls_ksp_type", "ksp type", "", NLS_KSP, NLS_KSP_TYPES, NLS_KSP[nlsP->ksp_type], &nlsP->ksp_type, 0);CHKERRQ(ierr);
966a7e14dcfSSatish 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);
967a7e14dcfSSatish 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);
968a7e14dcfSSatish 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);
969a7e14dcfSSatish 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);
970a7e14dcfSSatish Balay  ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval, 0);CHKERRQ(ierr);
971a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin, 0);CHKERRQ(ierr);
972a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax, 0);CHKERRQ(ierr);
973a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac, 0);CHKERRQ(ierr);
974a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin, 0);CHKERRQ(ierr);
975a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax, 0);CHKERRQ(ierr);
976a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac, 0);CHKERRQ(ierr);
977a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac, 0);CHKERRQ(ierr);
978a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac, 0);CHKERRQ(ierr);
979a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac, 0);CHKERRQ(ierr);
980a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1, 0);CHKERRQ(ierr);
981a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2, 0);CHKERRQ(ierr);
982a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3, 0);CHKERRQ(ierr);
983a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4, 0);CHKERRQ(ierr);
984a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1, 0);CHKERRQ(ierr);
985a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2, 0);CHKERRQ(ierr);
986a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3, 0);CHKERRQ(ierr);
987a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4, 0);CHKERRQ(ierr);
988a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5, 0);CHKERRQ(ierr);
989a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1, 0);CHKERRQ(ierr);
990a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2, 0);CHKERRQ(ierr);
991a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3, 0);CHKERRQ(ierr);
992a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4, 0);CHKERRQ(ierr);
993a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1, 0);CHKERRQ(ierr);
994a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2, 0);CHKERRQ(ierr);
995a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3, 0);CHKERRQ(ierr);
996a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4, 0);CHKERRQ(ierr);
997a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5, 0);CHKERRQ(ierr);
998a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i, 0);CHKERRQ(ierr);
999a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i, 0);CHKERRQ(ierr);
1000a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i, 0);CHKERRQ(ierr);
1001a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i, 0);CHKERRQ(ierr);
1002a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i, 0);CHKERRQ(ierr);
1003a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i, 0);CHKERRQ(ierr);
1004a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i, 0);CHKERRQ(ierr);
1005a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1, 0);CHKERRQ(ierr);
1006a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2, 0);CHKERRQ(ierr);
1007a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1, 0);CHKERRQ(ierr);
1008a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2, 0);CHKERRQ(ierr);
1009a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3, 0);CHKERRQ(ierr);
1010a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4, 0);CHKERRQ(ierr);
1011a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta, 0);CHKERRQ(ierr);
1012a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius, 0);CHKERRQ(ierr);
1013a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius, 0);CHKERRQ(ierr);
1014a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon, 0);CHKERRQ(ierr);
1015a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
1016a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1017a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1018a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1019a7e14dcfSSatish Balay }
1020a7e14dcfSSatish Balay 
1021a7e14dcfSSatish Balay 
1022a7e14dcfSSatish Balay /*------------------------------------------------------------*/
1023a7e14dcfSSatish Balay #undef __FUNCT__
1024a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NLS"
1025a7e14dcfSSatish Balay static PetscErrorCode TaoView_NLS(TaoSolver tao, PetscViewer viewer)
1026a7e14dcfSSatish Balay {
1027a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
1028a7e14dcfSSatish Balay   PetscInt       nrejects;
1029a7e14dcfSSatish Balay   PetscBool      isascii;
1030a7e14dcfSSatish Balay   PetscErrorCode ierr;
1031a7e14dcfSSatish Balay 
1032a7e14dcfSSatish Balay   PetscFunctionBegin;
1033a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1034a7e14dcfSSatish Balay   if (isascii) {
1035a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1036a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
1037a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr);
1038a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1039a7e14dcfSSatish Balay     }
1040a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
1041a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
1042a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
1043a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
1044a7e14dcfSSatish Balay 
1045a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
1046a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
1047a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
1048a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
1049a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
1050a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
1051a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
1052a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1053a7e14dcfSSatish Balay   }
1054a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1055a7e14dcfSSatish Balay }
1056a7e14dcfSSatish Balay 
1057a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
1058a7e14dcfSSatish Balay EXTERN_C_BEGIN
1059a7e14dcfSSatish Balay #undef __FUNCT__
1060a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NLS"
1061a7e14dcfSSatish Balay PetscErrorCode TaoCreate_NLS(TaoSolver tao)
1062a7e14dcfSSatish Balay {
1063a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
1064a7e14dcfSSatish Balay   const char     *morethuente_type = TAOLINESEARCH_MT;
1065a7e14dcfSSatish Balay   PetscErrorCode ierr;
1066a7e14dcfSSatish Balay 
1067a7e14dcfSSatish Balay   PetscFunctionBegin;
10683c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1069a7e14dcfSSatish Balay 
1070a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
1071a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
1072a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
1073a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1074a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
1075a7e14dcfSSatish Balay 
1076a7e14dcfSSatish Balay   tao->max_it = 50;
1077a7e14dcfSSatish Balay   tao->fatol = 1e-10;
1078a7e14dcfSSatish Balay   tao->frtol = 1e-10;
1079a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
1080a7e14dcfSSatish Balay   tao->trust0 = 100.0;
1081a7e14dcfSSatish Balay 
1082a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
1083a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
1084a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
1085a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
1086a7e14dcfSSatish Balay 
1087a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
1088a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
1089a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
1090a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
1091a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
1092a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
1093a7e14dcfSSatish Balay 
1094a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
1095a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
1096a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
1097a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
1098a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
1099a7e14dcfSSatish Balay 
1100a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
1101a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
1102a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
1103a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
1104a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
1105a7e14dcfSSatish Balay 
1106a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
1107a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
1108a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
1109a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
1110a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
1111a7e14dcfSSatish Balay 
1112a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
1113a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
1114a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
1115a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
1116a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
1117a7e14dcfSSatish Balay 
1118a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
1119a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
1120a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
1121a7e14dcfSSatish Balay 
1122a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
1123a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
1124a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
1125a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
1126a7e14dcfSSatish Balay 
1127a7e14dcfSSatish Balay   nlsP->theta = 0.05;
1128a7e14dcfSSatish Balay 
1129a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
1130a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
1131a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
1132a7e14dcfSSatish Balay 
1133a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
1134a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
1135a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
1136a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
1137a7e14dcfSSatish Balay 
1138a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
1139a7e14dcfSSatish Balay 
1140a7e14dcfSSatish Balay   /*  Remaining parameters */
1141a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
1142a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
1143a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
1144a7e14dcfSSatish Balay 
1145a7e14dcfSSatish Balay   nlsP->ksp_type        = NLS_KSP_STCG;
1146a7e14dcfSSatish Balay   nlsP->pc_type         = NLS_PC_BFGS;
1147a7e14dcfSSatish Balay   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1148a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1149a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
1150a7e14dcfSSatish Balay 
1151a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1152a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1153a7e14dcfSSatish Balay   ierr = TaoLineSearchUseTaoSolverRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1154a7e14dcfSSatish Balay 
1155a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
1156a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1157a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1158a7e14dcfSSatish Balay }
1159a7e14dcfSSatish Balay EXTERN_C_END
1160a7e14dcfSSatish Balay 
1161a7e14dcfSSatish Balay #undef __FUNCT__
1162a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
1163a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
1164a7e14dcfSSatish Balay {
1165a7e14dcfSSatish Balay   PetscErrorCode ierr;
1166a7e14dcfSSatish Balay   Mat            M;
116787f595a5SBarry Smith 
1168a7e14dcfSSatish Balay   PetscFunctionBegin;
1169a7e14dcfSSatish Balay   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1170a7e14dcfSSatish Balay   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
1171a7e14dcfSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
1172a7e14dcfSSatish Balay   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
1173a7e14dcfSSatish Balay   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
1174a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1175a7e14dcfSSatish Balay }
1176