xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision 2d9aa51bd4e54bdb7dc3042b9d4a676eba5efc58)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
3aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/nls/nls.h>
4a7e14dcfSSatish Balay 
5aaa7dc30SBarry Smith #include <petscksp.h>
6aaa7dc30SBarry Smith #include <petscpc.h>
7aaa7dc30SBarry Smith #include <petsc-private/kspimpl.h>
8aaa7dc30SBarry 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"
71441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao 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;
77e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
78e4cb33bbSBarry Smith   TaoConvergedReason           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   PetscInt                     stepType;
88a7e14dcfSSatish Balay   PetscInt                     iter = 0;
89a7e14dcfSSatish Balay   PetscInt                     bfgsUpdates = 0;
90a7e14dcfSSatish Balay   PetscInt                     n,N,kspits;
91a7e14dcfSSatish Balay   PetscInt                     needH;
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay   PetscInt                     i_max = 5;
94a7e14dcfSSatish Balay   PetscInt                     j_max = 1;
95a7e14dcfSSatish Balay   PetscInt                     i, j;
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay   PetscFunctionBegin;
98a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
99a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
100a7e14dcfSSatish Balay   }
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay   /* Initialized variables */
103a7e14dcfSSatish Balay   pert = nlsP->sval;
104a7e14dcfSSatish Balay 
105*2d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
106a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
107a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
108a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
109a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
110a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
111a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
112a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay   /* Modify the linear solver to a trust region method if desired */
115a7e14dcfSSatish Balay   switch(nlsP->ksp_type) {
116a7e14dcfSSatish Balay   case NLS_KSP_CG:
117a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr);
11872b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
119a7e14dcfSSatish Balay     break;
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   case NLS_KSP_NASH:
122a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
12372b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
124a7e14dcfSSatish Balay     break;
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay   case NLS_KSP_STCG:
127a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
12872b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
129a7e14dcfSSatish Balay     break;
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay   case NLS_KSP_GLTR:
132a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
13372b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
134a7e14dcfSSatish Balay     break;
135a7e14dcfSSatish Balay 
136a7e14dcfSSatish Balay   default:
137a7e14dcfSSatish Balay     /* Use the method set by the ksp_type */
138a7e14dcfSSatish Balay     break;
139a7e14dcfSSatish Balay   }
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
142a7e14dcfSSatish Balay    Will be reset during the first iteration */
143a7e14dcfSSatish Balay   if (NLS_KSP_NASH == nlsP->ksp_type) {
144a7e14dcfSSatish Balay     ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
145a7e14dcfSSatish Balay   } else if (NLS_KSP_STCG == nlsP->ksp_type) {
146a7e14dcfSSatish Balay     ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
147a7e14dcfSSatish Balay   } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
148a7e14dcfSSatish Balay     ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
149a7e14dcfSSatish Balay   }
150a7e14dcfSSatish Balay 
15187f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
152a7e14dcfSSatish Balay     tao->trust = tao->trust0;
153a7e14dcfSSatish Balay 
15487f595a5SBarry Smith     if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative");
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay     /* Modify the radius if it is too large or small */
157a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
158a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
159a7e14dcfSSatish Balay   }
160a7e14dcfSSatish Balay 
161a7e14dcfSSatish Balay   /* Get vectors we will need */
162a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
163a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
164a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
165a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
166a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr);
167a7e14dcfSSatish Balay   }
168a7e14dcfSSatish Balay 
169a7e14dcfSSatish Balay   /* Check convergence criteria */
170a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
171a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
17287f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
173a7e14dcfSSatish Balay   needH = 1;
174a7e14dcfSSatish Balay 
175a7e14dcfSSatish Balay   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
17687f595a5SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay   /* create vectors for the limited memory preconditioner */
17987f595a5SBarry Smith   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
180a7e14dcfSSatish Balay     if (!nlsP->Diag) {
181a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
182a7e14dcfSSatish Balay     }
183a7e14dcfSSatish Balay   }
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay   /* Modify the preconditioner to use the bfgs approximation */
186a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
187a7e14dcfSSatish Balay   switch(nlsP->pc_type) {
188a7e14dcfSSatish Balay   case NLS_PC_NONE:
189a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
190a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
191a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
192a7e14dcfSSatish Balay     }
193a7e14dcfSSatish Balay     break;
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay   case NLS_PC_AHESS:
196a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
197a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
198a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
199a7e14dcfSSatish Balay     }
200a7e14dcfSSatish Balay     ierr = PCJacobiSetUseAbs(pc);CHKERRQ(ierr);
201a7e14dcfSSatish Balay     break;
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay   case NLS_PC_BFGS:
204a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
205a7e14dcfSSatish Balay     if (pc->ops->setfromoptions) {
206a7e14dcfSSatish Balay       (*pc->ops->setfromoptions)(pc);
207a7e14dcfSSatish Balay     }
208a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
209a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
210a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
211a7e14dcfSSatish Balay     break;
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay   default:
214a7e14dcfSSatish Balay     /* Use the pc method set by pc_type */
215a7e14dcfSSatish Balay     break;
216a7e14dcfSSatish Balay   }
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
219a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
22087f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
221a7e14dcfSSatish Balay     switch(nlsP->init_type) {
222a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
223a7e14dcfSSatish Balay       /* Use the initial radius specified */
224a7e14dcfSSatish Balay       break;
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
227a7e14dcfSSatish Balay       /* Use the initial radius specified */
228a7e14dcfSSatish Balay       max_radius = 0.0;
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
231a7e14dcfSSatish Balay         fmin = f;
232a7e14dcfSSatish Balay         sigma = 0.0;
233a7e14dcfSSatish Balay 
234a7e14dcfSSatish Balay         if (needH) {
235ffad9901SBarry Smith           ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
236a7e14dcfSSatish Balay           needH = 0;
237a7e14dcfSSatish Balay         }
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
240a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
241a7e14dcfSSatish Balay           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
242a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
243a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
244a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
24587f595a5SBarry Smith           } else {
246a7e14dcfSSatish Balay             if (ftrial < fmin) {
247a7e14dcfSSatish Balay               fmin = ftrial;
248a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
249a7e14dcfSSatish Balay             }
250a7e14dcfSSatish Balay 
251a7e14dcfSSatish Balay             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
252a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
253a7e14dcfSSatish Balay 
254a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
255a7e14dcfSSatish Balay             actred = f - ftrial;
25687f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
257a7e14dcfSSatish Balay               kappa = 1.0;
25887f595a5SBarry Smith             } else {
259a7e14dcfSSatish Balay               kappa = actred / prered;
260a7e14dcfSSatish Balay             }
261a7e14dcfSSatish Balay 
262a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
263a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
264a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
265a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
266a7e14dcfSSatish Balay 
267a7e14dcfSSatish Balay             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
268a7e14dcfSSatish Balay               /* Great agreement */
269a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay               if (tau_max < 1.0) {
272a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
27387f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
274a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
27587f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
276a7e14dcfSSatish Balay                 tau = tau_1;
27787f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
278a7e14dcfSSatish Balay                 tau = tau_2;
27987f595a5SBarry Smith               } else {
280a7e14dcfSSatish Balay                 tau = tau_max;
281a7e14dcfSSatish Balay               }
28287f595a5SBarry Smith             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
283a7e14dcfSSatish Balay               /* Good agreement */
284a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
285a7e14dcfSSatish Balay 
286a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
287a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
28887f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
289a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
29087f595a5SBarry Smith               } else {
291a7e14dcfSSatish Balay                 tau = tau_max;
292a7e14dcfSSatish Balay               }
29387f595a5SBarry Smith             } else {
294a7e14dcfSSatish Balay               /* Not good agreement */
295a7e14dcfSSatish Balay               if (tau_min > 1.0) {
296a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
29787f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
298a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
29987f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
300a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
30187f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
302a7e14dcfSSatish Balay                 tau = tau_1;
30387f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
304a7e14dcfSSatish Balay                 tau = tau_2;
30587f595a5SBarry Smith               } else {
306a7e14dcfSSatish Balay                 tau = tau_max;
307a7e14dcfSSatish Balay               }
308a7e14dcfSSatish Balay             }
309a7e14dcfSSatish Balay           }
310a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
311a7e14dcfSSatish Balay         }
312a7e14dcfSSatish Balay 
313a7e14dcfSSatish Balay         if (fmin < f) {
314a7e14dcfSSatish Balay           f = fmin;
315a7e14dcfSSatish Balay           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
316a7e14dcfSSatish Balay           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
317a7e14dcfSSatish Balay 
318a7e14dcfSSatish Balay           ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
31987f595a5SBarry Smith           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
320a7e14dcfSSatish Balay           needH = 1;
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay           ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
32387f595a5SBarry Smith           if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
324a7e14dcfSSatish Balay         }
325a7e14dcfSSatish Balay       }
326a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
327a7e14dcfSSatish Balay 
328a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
329a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
330a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
331a7e14dcfSSatish Balay       break;
332a7e14dcfSSatish Balay 
333a7e14dcfSSatish Balay     default:
334a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
335a7e14dcfSSatish Balay       tao->trust = 0.0;
336a7e14dcfSSatish Balay       break;
337a7e14dcfSSatish Balay     }
338a7e14dcfSSatish Balay   }
339a7e14dcfSSatish Balay 
340a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
341a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
342a7e14dcfSSatish Balay      since the function value may have decreased */
343a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type) {
344a7e14dcfSSatish Balay     if (f != 0.0) {
345a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
34687f595a5SBarry Smith     } else {
347a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
348a7e14dcfSSatish Balay     }
349a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
350a7e14dcfSSatish Balay   }
351a7e14dcfSSatish Balay 
352a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
353a7e14dcfSSatish Balay   nlsP->newt = 0;
354a7e14dcfSSatish Balay   nlsP->bfgs = 0;
355a7e14dcfSSatish Balay   nlsP->sgrad = 0;
356a7e14dcfSSatish Balay   nlsP->grad = 0;
357a7e14dcfSSatish Balay 
358a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
359a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
360a7e14dcfSSatish Balay     ++iter;
361ae93cb3cSJason Sarich     tao->ksp_its=0;
362a7e14dcfSSatish Balay 
363a7e14dcfSSatish Balay     /* Compute the Hessian */
364a7e14dcfSSatish Balay     if (needH) {
365ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
366a7e14dcfSSatish Balay       needH = 0;
367a7e14dcfSSatish Balay     }
368a7e14dcfSSatish Balay 
36987f595a5SBarry Smith     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
370a7e14dcfSSatish Balay       /* Obtain diagonal for the bfgs preconditioner  */
371a7e14dcfSSatish Balay       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
372a7e14dcfSSatish Balay       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
373a7e14dcfSSatish Balay       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
374a7e14dcfSSatish Balay       ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
375a7e14dcfSSatish Balay     }
376a7e14dcfSSatish Balay 
377a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
378a7e14dcfSSatish Balay     if (pert > 0) {
379a7e14dcfSSatish Balay       ierr = MatShift(tao->hessian, pert);
380a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
381a7e14dcfSSatish Balay         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
382a7e14dcfSSatish Balay       }
383a7e14dcfSSatish Balay     }
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type) {
386a7e14dcfSSatish Balay       if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) {
387a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner  */
388a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
389a7e14dcfSSatish Balay         ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
390a7e14dcfSSatish Balay         ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
391a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
392a7e14dcfSSatish Balay       }
393a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
394a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
395a7e14dcfSSatish Balay       ++bfgsUpdates;
396a7e14dcfSSatish Balay     }
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
39923ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
40087f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type ||  NLS_KSP_GLTR == nlsP->ksp_type) {
401a7e14dcfSSatish Balay 
402a7e14dcfSSatish Balay       if (NLS_KSP_NASH == nlsP->ksp_type) {
403a7e14dcfSSatish Balay         ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
404a7e14dcfSSatish Balay       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
405a7e14dcfSSatish Balay          ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
406a7e14dcfSSatish Balay       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
407a7e14dcfSSatish Balay         ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
408a7e14dcfSSatish Balay       }
409a7e14dcfSSatish Balay 
410a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
411a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
412a7e14dcfSSatish Balay       tao->ksp_its+=kspits;
413ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
414a7e14dcfSSatish Balay 
415a7e14dcfSSatish Balay       if (NLS_KSP_NASH == nlsP->ksp_type) {
416a7e14dcfSSatish Balay         ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
417a7e14dcfSSatish Balay       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
418a7e14dcfSSatish Balay          ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
419a7e14dcfSSatish Balay       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
420a7e14dcfSSatish Balay         ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
421a7e14dcfSSatish Balay       }
422a7e14dcfSSatish Balay 
423a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
424a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
425a7e14dcfSSatish Balay         if (norm_d > 0.0) {
426a7e14dcfSSatish Balay           tao->trust = norm_d;
427a7e14dcfSSatish Balay 
428a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
429a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
430a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
43187f595a5SBarry Smith         } else {
432a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
433a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
434a7e14dcfSSatish Balay           tao->trust = tao->trust0;
435a7e14dcfSSatish Balay 
436a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
437a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
438a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
439a7e14dcfSSatish Balay 
440a7e14dcfSSatish Balay           if (NLS_KSP_NASH == nlsP->ksp_type) {
441a7e14dcfSSatish Balay             ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
442a7e14dcfSSatish Balay           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
443a7e14dcfSSatish Balay             ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
444a7e14dcfSSatish Balay           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
445a7e14dcfSSatish Balay             ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
446a7e14dcfSSatish Balay           }
447a7e14dcfSSatish Balay 
448a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
449a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
450a7e14dcfSSatish Balay           tao->ksp_its+=kspits;
451ae93cb3cSJason Sarich           tao->ksp_tot_its+=kspits;
452a7e14dcfSSatish Balay           if (NLS_KSP_NASH == nlsP->ksp_type) {
453a7e14dcfSSatish Balay             ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
454a7e14dcfSSatish Balay           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
455a7e14dcfSSatish Balay             ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
456a7e14dcfSSatish Balay           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
457a7e14dcfSSatish Balay             ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
458a7e14dcfSSatish Balay           }
459a7e14dcfSSatish Balay 
46087f595a5SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
461a7e14dcfSSatish Balay         }
462a7e14dcfSSatish Balay       }
46387f595a5SBarry Smith     } else {
464a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
465a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
466a7e14dcfSSatish Balay       tao->ksp_its += kspits;
467ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
468a7e14dcfSSatish Balay     }
469a7e14dcfSSatish Balay     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
470a7e14dcfSSatish Balay     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
47187f595a5SBarry Smith     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
472a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
473a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
474a7e14dcfSSatish Balay 
475a7e14dcfSSatish Balay       if (f != 0.0) {
476a7e14dcfSSatish Balay         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
47787f595a5SBarry Smith       } else {
478a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
479a7e14dcfSSatish Balay       }
480a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
481a7e14dcfSSatish Balay       ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
482a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
483a7e14dcfSSatish Balay       bfgsUpdates = 1;
484a7e14dcfSSatish Balay     }
485a7e14dcfSSatish Balay 
486a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
487a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
48887f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
489a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
49087f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
491a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
49287f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
493a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
49487f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
495a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
49687f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
497a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
49887f595a5SBarry Smith     } else {
499a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
500a7e14dcfSSatish Balay     }
501a7e14dcfSSatish Balay 
502a7e14dcfSSatish Balay     /* Check for success (descent direction) */
503a7e14dcfSSatish Balay     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
504a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
505a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
506a7e14dcfSSatish Balay          Update the perturbation for next time */
507a7e14dcfSSatish Balay       if (pert <= 0.0) {
508a7e14dcfSSatish Balay         /* Initialize the perturbation */
509a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
510a7e14dcfSSatish Balay         if (NLS_KSP_GLTR == nlsP->ksp_type) {
511a7e14dcfSSatish Balay           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
512a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
513a7e14dcfSSatish Balay         }
51487f595a5SBarry Smith       } else {
515a7e14dcfSSatish Balay         /* Increase the perturbation */
516a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
517a7e14dcfSSatish Balay       }
518a7e14dcfSSatish Balay 
519a7e14dcfSSatish Balay       if (NLS_PC_BFGS != nlsP->pc_type) {
520a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
521a7e14dcfSSatish Balay            Must use gradient direction in this case */
522a7e14dcfSSatish Balay         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
523a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
524a7e14dcfSSatish Balay         ++nlsP->grad;
525a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
52687f595a5SBarry Smith       } else {
527a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
528a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
529a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
530a7e14dcfSSatish Balay 
531a7e14dcfSSatish Balay         /* Check for success (descent direction) */
532a7e14dcfSSatish Balay         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
533a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
534a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
535a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
536a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
537a7e14dcfSSatish Balay              which is guaranteed to be descent */
538a7e14dcfSSatish Balay 
539a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
540a7e14dcfSSatish Balay 
541a7e14dcfSSatish Balay           if (f != 0.0) {
542a7e14dcfSSatish Balay             delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
54387f595a5SBarry Smith           } else {
544a7e14dcfSSatish Balay             delta = 2.0 / (gnorm*gnorm);
545a7e14dcfSSatish Balay           }
546a7e14dcfSSatish Balay           ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
547a7e14dcfSSatish Balay           ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
548a7e14dcfSSatish Balay           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
549a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
550a7e14dcfSSatish Balay           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
551a7e14dcfSSatish Balay 
552a7e14dcfSSatish Balay           bfgsUpdates = 1;
553a7e14dcfSSatish Balay           ++nlsP->sgrad;
554a7e14dcfSSatish Balay           stepType = NLS_SCALED_GRADIENT;
55587f595a5SBarry Smith         } else {
556a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
557a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
558a7e14dcfSSatish Balay             ++nlsP->sgrad;
559a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
56087f595a5SBarry Smith           } else {
561a7e14dcfSSatish Balay             ++nlsP->bfgs;
562a7e14dcfSSatish Balay             stepType = NLS_BFGS;
563a7e14dcfSSatish Balay           }
564a7e14dcfSSatish Balay         }
565a7e14dcfSSatish Balay       }
56687f595a5SBarry Smith     } else {
567a7e14dcfSSatish Balay       /* Computed Newton step is descent */
568a7e14dcfSSatish Balay       switch (ksp_reason) {
569a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
570a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
571a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
572a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
573a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
574a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
575a7e14dcfSSatish Balay         if (pert <= 0.0) {
576a7e14dcfSSatish Balay           /* Initialize the perturbation */
577a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
578a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
579a7e14dcfSSatish Balay             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
580a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
581a7e14dcfSSatish Balay           }
58287f595a5SBarry Smith         } else {
583a7e14dcfSSatish Balay           /* Increase the perturbation */
584a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
585a7e14dcfSSatish Balay         }
586a7e14dcfSSatish Balay         break;
587a7e14dcfSSatish Balay 
588a7e14dcfSSatish Balay       default:
589a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
590a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
591a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
592a7e14dcfSSatish Balay           pert = 0.0;
593a7e14dcfSSatish Balay         }
594a7e14dcfSSatish Balay         break;
595a7e14dcfSSatish Balay       }
596a7e14dcfSSatish Balay 
597a7e14dcfSSatish Balay       ++nlsP->newt;
598a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
599a7e14dcfSSatish Balay     }
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay     /* Perform the linesearch */
602a7e14dcfSSatish Balay     fold = f;
603a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
604a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
605a7e14dcfSSatish Balay 
606a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
607a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
608a7e14dcfSSatish Balay 
60987f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
610a7e14dcfSSatish Balay       f = fold;
611a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
612a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
613a7e14dcfSSatish Balay 
614a7e14dcfSSatish Balay       switch(stepType) {
615a7e14dcfSSatish Balay       case NLS_NEWTON:
616a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
617a7e14dcfSSatish Balay            Update the perturbation for next time */
618a7e14dcfSSatish Balay         if (pert <= 0.0) {
619a7e14dcfSSatish Balay           /* Initialize the perturbation */
620a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
621a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
622a7e14dcfSSatish Balay             ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
623a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
624a7e14dcfSSatish Balay           }
62587f595a5SBarry Smith         } else {
626a7e14dcfSSatish Balay           /* Increase the perturbation */
627a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
628a7e14dcfSSatish Balay         }
629a7e14dcfSSatish Balay 
630a7e14dcfSSatish Balay         if (NLS_PC_BFGS != nlsP->pc_type) {
631a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
632a7e14dcfSSatish Balay              Must use gradient direction in this case */
633a7e14dcfSSatish Balay           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
634a7e14dcfSSatish Balay           ++nlsP->grad;
635a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
63687f595a5SBarry Smith         } else {
637a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
638a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
639a7e14dcfSSatish Balay           /* Check for success (descent direction) */
640a7e14dcfSSatish Balay           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
641a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
642a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
643a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
644a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
645a7e14dcfSSatish Balay 
646a7e14dcfSSatish Balay             if (f != 0.0) {
647a7e14dcfSSatish Balay               delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
64887f595a5SBarry Smith             } else {
649a7e14dcfSSatish Balay               delta = 2.0 / (gnorm*gnorm);
650a7e14dcfSSatish Balay             }
651a7e14dcfSSatish Balay             ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
652a7e14dcfSSatish Balay             ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
653a7e14dcfSSatish Balay             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
654a7e14dcfSSatish Balay             ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
655a7e14dcfSSatish Balay 
656a7e14dcfSSatish Balay             bfgsUpdates = 1;
657a7e14dcfSSatish Balay             ++nlsP->sgrad;
658a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
65987f595a5SBarry Smith           } else {
660a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
661a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
662a7e14dcfSSatish Balay               ++nlsP->sgrad;
663a7e14dcfSSatish Balay               stepType = NLS_SCALED_GRADIENT;
66487f595a5SBarry Smith             } else {
665a7e14dcfSSatish Balay               ++nlsP->bfgs;
666a7e14dcfSSatish Balay               stepType = NLS_BFGS;
667a7e14dcfSSatish Balay             }
668a7e14dcfSSatish Balay           }
669a7e14dcfSSatish Balay         }
670a7e14dcfSSatish Balay         break;
671a7e14dcfSSatish Balay 
672a7e14dcfSSatish Balay       case NLS_BFGS:
673a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
674a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
675a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
676a7e14dcfSSatish Balay 
677a7e14dcfSSatish Balay         if (f != 0.0) {
678a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
67987f595a5SBarry Smith         } else {
680a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
681a7e14dcfSSatish Balay         }
682a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
683a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
684a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
685a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
686a7e14dcfSSatish Balay 
687a7e14dcfSSatish Balay         bfgsUpdates = 1;
688a7e14dcfSSatish Balay         ++nlsP->sgrad;
689a7e14dcfSSatish Balay         stepType = NLS_SCALED_GRADIENT;
690a7e14dcfSSatish Balay         break;
691a7e14dcfSSatish Balay 
692a7e14dcfSSatish Balay       case NLS_SCALED_GRADIENT:
693a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
694a7e14dcfSSatish Balay            The scaled gradient step did not produce a new iterate;
695a7e14dcfSSatish Balay            attemp to use the gradient direction.
696a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
697a7e14dcfSSatish Balay 
698a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr);
699a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr);
700a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
701a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
702a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
703a7e14dcfSSatish Balay 
704a7e14dcfSSatish Balay         bfgsUpdates = 1;
705a7e14dcfSSatish Balay         ++nlsP->grad;
706a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
707a7e14dcfSSatish Balay         break;
708a7e14dcfSSatish Balay       }
709a7e14dcfSSatish Balay       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
710a7e14dcfSSatish Balay 
711a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
712a7e14dcfSSatish Balay       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
713a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
714a7e14dcfSSatish Balay     }
715a7e14dcfSSatish Balay 
71687f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
717a7e14dcfSSatish Balay       /* Failed to find an improving point */
718a7e14dcfSSatish Balay       f = fold;
719a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
720a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
721a7e14dcfSSatish Balay       step = 0.0;
722a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
723a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
724a7e14dcfSSatish Balay       break;
725a7e14dcfSSatish Balay     }
726a7e14dcfSSatish Balay 
727a7e14dcfSSatish Balay     /* Update trust region radius */
72887f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
729a7e14dcfSSatish Balay       switch(nlsP->update_type) {
730a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
731a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
732a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
733a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
734a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
73587f595a5SBarry Smith           } else if (step < nlsP->nu2) {
736a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
737a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
73887f595a5SBarry Smith           } else if (step < nlsP->nu3) {
739a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
740a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
741a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
74287f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
743a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
744a7e14dcfSSatish Balay             }
74587f595a5SBarry Smith           } else if (step < nlsP->nu4) {
746a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
747a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
74887f595a5SBarry Smith           } else {
749a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
750a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
751a7e14dcfSSatish Balay           }
75287f595a5SBarry Smith         } else {
753a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
754a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
755a7e14dcfSSatish Balay         }
756a7e14dcfSSatish Balay         break;
757a7e14dcfSSatish Balay 
758a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
759a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
760a7e14dcfSSatish Balay           /*  Get predicted reduction */
761a7e14dcfSSatish Balay 
762a7e14dcfSSatish Balay           if (NLS_KSP_STCG == nlsP->ksp_type) {
763a7e14dcfSSatish Balay             ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);
764a7e14dcfSSatish Balay           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
765a7e14dcfSSatish Balay             ierr = KSPNASHGetObjFcn(tao->ksp,&prered);
766a7e14dcfSSatish Balay           } else {
767a7e14dcfSSatish Balay             ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
768a7e14dcfSSatish Balay           }
769a7e14dcfSSatish Balay 
770a7e14dcfSSatish Balay           if (prered >= 0.0) {
771a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
772a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
773a7e14dcfSSatish Balay             /*  be rejected! */
774a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
77587f595a5SBarry Smith           } else {
776a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
777a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
77887f595a5SBarry Smith             } else {
779a7e14dcfSSatish Balay               /*  Compute and actual reduction */
780a7e14dcfSSatish Balay               actred = fold - f_full;
781a7e14dcfSSatish Balay               prered = -prered;
782a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
783a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
784a7e14dcfSSatish Balay                 kappa = 1.0;
78587f595a5SBarry Smith               } else {
786a7e14dcfSSatish Balay                 kappa = actred / prered;
787a7e14dcfSSatish Balay               }
788a7e14dcfSSatish Balay 
789a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
790a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
791a7e14dcfSSatish Balay                 /*  Very bad step */
792a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
79387f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
794a7e14dcfSSatish Balay                 /*  Marginal bad step */
795a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
79687f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
797a7e14dcfSSatish Balay                 /*  Reasonable step */
798a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
799a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
80087f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
801a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
802a7e14dcfSSatish Balay                 }
80387f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
804a7e14dcfSSatish Balay                 /*  Good step */
805a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
80687f595a5SBarry Smith               } else {
807a7e14dcfSSatish Balay                 /*  Very good step */
808a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
809a7e14dcfSSatish Balay               }
810a7e14dcfSSatish Balay             }
811a7e14dcfSSatish Balay           }
81287f595a5SBarry Smith         } else {
813a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
814a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
815a7e14dcfSSatish Balay         }
816a7e14dcfSSatish Balay         break;
817a7e14dcfSSatish Balay 
818a7e14dcfSSatish Balay       default:
819a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
820a7e14dcfSSatish Balay 
821a7e14dcfSSatish Balay           if (NLS_KSP_STCG == nlsP->ksp_type) {
822a7e14dcfSSatish Balay               ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);
823a7e14dcfSSatish Balay           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
824a7e14dcfSSatish Balay               ierr = KSPNASHGetObjFcn(tao->ksp,&prered);
825a7e14dcfSSatish Balay           } else {
826a7e14dcfSSatish Balay               ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
827a7e14dcfSSatish Balay           }
828a7e14dcfSSatish Balay           if (prered >= 0.0) {
829a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
830a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
831a7e14dcfSSatish Balay             /*  be rejected! */
832a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
83387f595a5SBarry Smith           } else {
834a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
835a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
83687f595a5SBarry Smith             } else {
837a7e14dcfSSatish Balay               actred = fold - f_full;
838a7e14dcfSSatish Balay               prered = -prered;
83987f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
840a7e14dcfSSatish Balay                 kappa = 1.0;
84187f595a5SBarry Smith               } else {
842a7e14dcfSSatish Balay                 kappa = actred / prered;
843a7e14dcfSSatish Balay               }
844a7e14dcfSSatish Balay 
845a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
846a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
847a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
848a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
849a7e14dcfSSatish Balay 
850a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
851a7e14dcfSSatish Balay                 /*  Great agreement */
852a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
853a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
85487f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
855a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
85687f595a5SBarry Smith                 } else {
857a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
858a7e14dcfSSatish Balay                 }
85987f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
860a7e14dcfSSatish Balay                 /*  Good agreement */
861a7e14dcfSSatish Balay 
862a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
863a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
86487f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
865a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
86687f595a5SBarry Smith                 } else if (tau_max < 1.0) {
867a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
86887f595a5SBarry Smith                 } else {
869a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
870a7e14dcfSSatish Balay                 }
87187f595a5SBarry Smith               } else {
872a7e14dcfSSatish Balay                 /*  Not good agreement */
873a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
874a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
87587f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
876a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
87787f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
878a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
87987f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
880a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
88187f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
882a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
88387f595a5SBarry Smith                 } else {
884a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
885a7e14dcfSSatish Balay                 }
886a7e14dcfSSatish Balay               }
887a7e14dcfSSatish Balay             }
888a7e14dcfSSatish Balay           }
88987f595a5SBarry Smith         } else {
890a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
891a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
892a7e14dcfSSatish Balay         }
893a7e14dcfSSatish Balay         break;
894a7e14dcfSSatish Balay       }
895a7e14dcfSSatish Balay 
896a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
897a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
898a7e14dcfSSatish Balay     }
899a7e14dcfSSatish Balay 
900a7e14dcfSSatish Balay     /*  Check for termination */
901a7e14dcfSSatish Balay     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
90287f595a5SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
903a7e14dcfSSatish Balay     needH = 1;
904a7e14dcfSSatish Balay     ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
905a7e14dcfSSatish Balay   }
906a7e14dcfSSatish Balay   PetscFunctionReturn(0);
907a7e14dcfSSatish Balay }
908a7e14dcfSSatish Balay 
909a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
910a7e14dcfSSatish Balay #undef __FUNCT__
911a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NLS"
912441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao)
913a7e14dcfSSatish Balay {
914a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
915a7e14dcfSSatish Balay   PetscErrorCode ierr;
916a7e14dcfSSatish Balay 
917a7e14dcfSSatish Balay   PetscFunctionBegin;
918a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
919a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
920a7e14dcfSSatish Balay   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
921a7e14dcfSSatish Balay   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
922a7e14dcfSSatish Balay   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
923a7e14dcfSSatish Balay   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
924a7e14dcfSSatish Balay   nlsP->Diag = 0;
925a7e14dcfSSatish Balay   nlsP->M = 0;
926a7e14dcfSSatish Balay   PetscFunctionReturn(0);
927a7e14dcfSSatish Balay }
928a7e14dcfSSatish Balay 
929a7e14dcfSSatish Balay /*------------------------------------------------------------*/
930a7e14dcfSSatish Balay #undef __FUNCT__
931a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NLS"
932441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao)
933a7e14dcfSSatish Balay {
934a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
935a7e14dcfSSatish Balay   PetscErrorCode ierr;
936a7e14dcfSSatish Balay 
937a7e14dcfSSatish Balay   PetscFunctionBegin;
938a7e14dcfSSatish Balay   if (tao->setupcalled) {
939a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
940a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
941a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
942a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
943a7e14dcfSSatish Balay   }
944a7e14dcfSSatish Balay   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
945a7e14dcfSSatish Balay   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
946a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
947a7e14dcfSSatish Balay   PetscFunctionReturn(0);
948a7e14dcfSSatish Balay }
949a7e14dcfSSatish Balay 
950a7e14dcfSSatish Balay /*------------------------------------------------------------*/
951a7e14dcfSSatish Balay #undef __FUNCT__
952a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NLS"
953441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(Tao tao)
954a7e14dcfSSatish Balay {
955a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
956a7e14dcfSSatish Balay   PetscErrorCode ierr;
957a7e14dcfSSatish Balay 
958a7e14dcfSSatish Balay   PetscFunctionBegin;
959a7e14dcfSSatish Balay   ierr = PetscOptionsHead("Newton line search method for unconstrained optimization");CHKERRQ(ierr);
960a7e14dcfSSatish 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);
961a7e14dcfSSatish 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);
962a7e14dcfSSatish 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);
963a7e14dcfSSatish 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);
964a7e14dcfSSatish 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);
965a7e14dcfSSatish Balay  ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval, 0);CHKERRQ(ierr);
966a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin, 0);CHKERRQ(ierr);
967a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax, 0);CHKERRQ(ierr);
968a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac, 0);CHKERRQ(ierr);
969a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin, 0);CHKERRQ(ierr);
970a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax, 0);CHKERRQ(ierr);
971a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac, 0);CHKERRQ(ierr);
972a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac, 0);CHKERRQ(ierr);
973a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac, 0);CHKERRQ(ierr);
974a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac, 0);CHKERRQ(ierr);
975a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1, 0);CHKERRQ(ierr);
976a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2, 0);CHKERRQ(ierr);
977a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3, 0);CHKERRQ(ierr);
978a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4, 0);CHKERRQ(ierr);
979a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1, 0);CHKERRQ(ierr);
980a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2, 0);CHKERRQ(ierr);
981a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3, 0);CHKERRQ(ierr);
982a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4, 0);CHKERRQ(ierr);
983a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5, 0);CHKERRQ(ierr);
984a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1, 0);CHKERRQ(ierr);
985a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2, 0);CHKERRQ(ierr);
986a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3, 0);CHKERRQ(ierr);
987a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4, 0);CHKERRQ(ierr);
988a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1, 0);CHKERRQ(ierr);
989a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2, 0);CHKERRQ(ierr);
990a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3, 0);CHKERRQ(ierr);
991a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4, 0);CHKERRQ(ierr);
992a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5, 0);CHKERRQ(ierr);
993a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i, 0);CHKERRQ(ierr);
994a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i, 0);CHKERRQ(ierr);
995a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i, 0);CHKERRQ(ierr);
996a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i, 0);CHKERRQ(ierr);
997a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i, 0);CHKERRQ(ierr);
998a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i, 0);CHKERRQ(ierr);
999a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i, 0);CHKERRQ(ierr);
1000a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1, 0);CHKERRQ(ierr);
1001a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2, 0);CHKERRQ(ierr);
1002a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1, 0);CHKERRQ(ierr);
1003a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2, 0);CHKERRQ(ierr);
1004a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3, 0);CHKERRQ(ierr);
1005a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4, 0);CHKERRQ(ierr);
1006a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta, 0);CHKERRQ(ierr);
1007a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius, 0);CHKERRQ(ierr);
1008a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius, 0);CHKERRQ(ierr);
1009a7e14dcfSSatish Balay   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon, 0);CHKERRQ(ierr);
1010a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
1011a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1012a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1013a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1014a7e14dcfSSatish Balay }
1015a7e14dcfSSatish Balay 
1016a7e14dcfSSatish Balay 
1017a7e14dcfSSatish Balay /*------------------------------------------------------------*/
1018a7e14dcfSSatish Balay #undef __FUNCT__
1019a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NLS"
1020441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
1021a7e14dcfSSatish Balay {
1022a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
1023a7e14dcfSSatish Balay   PetscInt       nrejects;
1024a7e14dcfSSatish Balay   PetscBool      isascii;
1025a7e14dcfSSatish Balay   PetscErrorCode ierr;
1026a7e14dcfSSatish Balay 
1027a7e14dcfSSatish Balay   PetscFunctionBegin;
1028a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1029a7e14dcfSSatish Balay   if (isascii) {
1030a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1031a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
1032a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr);
1033a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1034a7e14dcfSSatish Balay     }
1035a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
1036a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
1037a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
1038a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
1039a7e14dcfSSatish Balay 
1040a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
1041a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
1042a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
1043a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
1044a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
1045a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
1046a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
1047a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1048a7e14dcfSSatish Balay   }
1049a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1050a7e14dcfSSatish Balay }
1051a7e14dcfSSatish Balay 
1052a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
10534aa34175SJason Sarich /*MC
10544aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
10554aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
10564aa34175SJason Sarich   system of equations to obtain the step diretion dk:
10574aa34175SJason Sarich               Hk dk = -gk
10584aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
10594aa34175SJason Sarich   solve
10604aa34175SJason Sarich               min_t f(xk + t d_k)
10614aa34175SJason Sarich 
10624aa34175SJason Sarich     Options Database Keys:
10634aa34175SJason Sarich + -tao_nls_ksp_type - "cg","nash","stcg","gltr","petsc"
10644aa34175SJason Sarich . -tao_nls_pc_type - "none","ahess","bfgs","petsc"
10654aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs"
10664aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation"
10674aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
10684aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
10694aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
10704aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
10714aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
10724aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
10734aa34175SJason Sarich . -tao_nls_pgfac - growth factor
10744aa34175SJason Sarich . -tao_nls_psfac - shrink factor
10754aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
10764aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
10774aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
10784aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
10794aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
10804aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
10814aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
10824aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
10834aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
10844aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
10854aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
10864aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
10874aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
10884aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
10894aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
10904aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
10914aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
10924aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
10934aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
10944aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
10954aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
10964aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
10974aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
10984aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
10991522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
11001522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
11011522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
11021522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
11031522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
11041522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
11051522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
11061eb8069cSJason Sarich 
11071eb8069cSJason Sarich   Level: beginner
11081522df2eSJason Sarich M*/
11094aa34175SJason Sarich 
1110a7e14dcfSSatish Balay #undef __FUNCT__
1111a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NLS"
1112728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
1113a7e14dcfSSatish Balay {
1114a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
11158caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
1116a7e14dcfSSatish Balay   PetscErrorCode ierr;
1117a7e14dcfSSatish Balay 
1118a7e14dcfSSatish Balay   PetscFunctionBegin;
11193c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1120a7e14dcfSSatish Balay 
1121a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
1122a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
1123a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
1124a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1125a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
1126a7e14dcfSSatish Balay 
1127a7e14dcfSSatish Balay   tao->max_it = 50;
11286f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
11296f4723b1SBarry Smith   tao->fatol = 1e-5;
11306f4723b1SBarry Smith   tao->frtol = 1e-5;
11316f4723b1SBarry Smith #else
1132a7e14dcfSSatish Balay   tao->fatol = 1e-10;
1133a7e14dcfSSatish Balay   tao->frtol = 1e-10;
11346f4723b1SBarry Smith #endif
1135a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
1136a7e14dcfSSatish Balay   tao->trust0 = 100.0;
1137a7e14dcfSSatish Balay 
1138a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
1139a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
1140a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
1141a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
1142a7e14dcfSSatish Balay 
1143a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
1144a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
1145a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
1146a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
1147a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
1148a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
1149a7e14dcfSSatish Balay 
1150a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
1151a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
1152a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
1153a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
1154a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
1155a7e14dcfSSatish Balay 
1156a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
1157a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
1158a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
1159a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
1160a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
1161a7e14dcfSSatish Balay 
1162a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
1163a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
1164a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
1165a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
1166a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
1167a7e14dcfSSatish Balay 
1168a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
1169a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
1170a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
1171a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
1172a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
1173a7e14dcfSSatish Balay 
1174a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
1175a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
1176a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
1177a7e14dcfSSatish Balay 
1178a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
1179a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
1180a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
1181a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
1182a7e14dcfSSatish Balay 
1183a7e14dcfSSatish Balay   nlsP->theta = 0.05;
1184a7e14dcfSSatish Balay 
1185a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
1186a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
1187a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
1188a7e14dcfSSatish Balay 
1189a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
1190a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
1191a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
1192a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
1193a7e14dcfSSatish Balay 
1194a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
1195a7e14dcfSSatish Balay 
1196a7e14dcfSSatish Balay   /*  Remaining parameters */
1197a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
1198a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
1199a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
1200a7e14dcfSSatish Balay 
1201a7e14dcfSSatish Balay   nlsP->ksp_type        = NLS_KSP_STCG;
1202a7e14dcfSSatish Balay   nlsP->pc_type         = NLS_PC_BFGS;
1203a7e14dcfSSatish Balay   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1204a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1205a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
1206a7e14dcfSSatish Balay 
1207a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1208a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1209441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1210a7e14dcfSSatish Balay 
1211a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
1212a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1213a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1214a7e14dcfSSatish Balay }
1215a7e14dcfSSatish Balay 
1216a7e14dcfSSatish Balay #undef __FUNCT__
1217a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
1218a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
1219a7e14dcfSSatish Balay {
1220a7e14dcfSSatish Balay   PetscErrorCode ierr;
1221a7e14dcfSSatish Balay   Mat            M;
122287f595a5SBarry Smith 
1223a7e14dcfSSatish Balay   PetscFunctionBegin;
1224a7e14dcfSSatish Balay   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1225a7e14dcfSSatish Balay   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
1226a7e14dcfSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
1227a7e14dcfSSatish Balay   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
1228a7e14dcfSSatish Balay   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
1229a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1230a7e14dcfSSatish Balay }
1231