xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision a9603a14f409ff31cc951cd8cacb5ec045bdcdca)
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>
7af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
8af0996ceSBarry 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                     bfgsUpdates = 0;
89a7e14dcfSSatish Balay   PetscInt                     n,N,kspits;
90a7e14dcfSSatish Balay   PetscInt                     needH;
91a7e14dcfSSatish Balay 
92a7e14dcfSSatish Balay   PetscInt                     i_max = 5;
93a7e14dcfSSatish Balay   PetscInt                     j_max = 1;
94a7e14dcfSSatish Balay   PetscInt                     i, j;
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay   PetscFunctionBegin;
97a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
98a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
99a7e14dcfSSatish Balay   }
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay   /* Initialized variables */
102a7e14dcfSSatish Balay   pert = nlsP->sval;
103a7e14dcfSSatish Balay 
1042d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
105a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
106a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
107a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
108a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
109a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
110a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
111a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
112a7e14dcfSSatish Balay 
113a7e14dcfSSatish Balay   /* Modify the linear solver to a trust region method if desired */
114a7e14dcfSSatish Balay   switch(nlsP->ksp_type) {
115a7e14dcfSSatish Balay   case NLS_KSP_CG:
116a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr);
11772b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
118a7e14dcfSSatish Balay     break;
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay   case NLS_KSP_NASH:
121a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
12272b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
123a7e14dcfSSatish Balay     break;
124a7e14dcfSSatish Balay 
125a7e14dcfSSatish Balay   case NLS_KSP_STCG:
126a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
12772b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
128a7e14dcfSSatish Balay     break;
129a7e14dcfSSatish Balay 
130a7e14dcfSSatish Balay   case NLS_KSP_GLTR:
131a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
13272b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
133a7e14dcfSSatish Balay     break;
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay   default:
136a7e14dcfSSatish Balay     /* Use the method set by the ksp_type */
137a7e14dcfSSatish Balay     break;
138a7e14dcfSSatish Balay   }
139a7e14dcfSSatish Balay 
140a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
141a7e14dcfSSatish Balay    Will be reset during the first iteration */
142a7e14dcfSSatish Balay   if (NLS_KSP_NASH == nlsP->ksp_type) {
143a7e14dcfSSatish Balay     ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
144a7e14dcfSSatish Balay   } else if (NLS_KSP_STCG == nlsP->ksp_type) {
145a7e14dcfSSatish Balay     ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
146a7e14dcfSSatish Balay   } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
147a7e14dcfSSatish Balay     ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
148a7e14dcfSSatish Balay   }
149a7e14dcfSSatish Balay 
15087f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
151a7e14dcfSSatish Balay     tao->trust = tao->trust0;
152a7e14dcfSSatish Balay 
15387f595a5SBarry Smith     if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative");
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay     /* Modify the radius if it is too large or small */
156a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
157a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
158a7e14dcfSSatish Balay   }
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay   /* Get vectors we will need */
161a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
162a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
163a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
164a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
165a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr);
166a7e14dcfSSatish Balay   }
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   /* Check convergence criteria */
169a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
170*a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
17187f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
172a7e14dcfSSatish Balay   needH = 1;
173a7e14dcfSSatish Balay 
1748931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
17587f595a5SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
176a7e14dcfSSatish Balay 
177a7e14dcfSSatish Balay   /* create vectors for the limited memory preconditioner */
17887f595a5SBarry Smith   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
179a7e14dcfSSatish Balay     if (!nlsP->Diag) {
180a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
181a7e14dcfSSatish Balay     }
182a7e14dcfSSatish Balay   }
183a7e14dcfSSatish Balay 
184a7e14dcfSSatish Balay   /* Modify the preconditioner to use the bfgs approximation */
185a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
186a7e14dcfSSatish Balay   switch(nlsP->pc_type) {
187a7e14dcfSSatish Balay   case NLS_PC_NONE:
188a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1891a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
190a7e14dcfSSatish Balay     break;
191a7e14dcfSSatish Balay 
192a7e14dcfSSatish Balay   case NLS_PC_AHESS:
193a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1941a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
195baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
196a7e14dcfSSatish Balay     break;
197a7e14dcfSSatish Balay 
198a7e14dcfSSatish Balay   case NLS_PC_BFGS:
199a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
2001a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
201a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
202a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
203a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
204a7e14dcfSSatish Balay     break;
205a7e14dcfSSatish Balay 
206a7e14dcfSSatish Balay   default:
207a7e14dcfSSatish Balay     /* Use the pc method set by pc_type */
208a7e14dcfSSatish Balay     break;
209a7e14dcfSSatish Balay   }
210a7e14dcfSSatish Balay 
211a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
212a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
21387f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
214a7e14dcfSSatish Balay     switch(nlsP->init_type) {
215a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
216a7e14dcfSSatish Balay       /* Use the initial radius specified */
217a7e14dcfSSatish Balay       break;
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
220a7e14dcfSSatish Balay       /* Use the initial radius specified */
221a7e14dcfSSatish Balay       max_radius = 0.0;
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
224a7e14dcfSSatish Balay         fmin = f;
225a7e14dcfSSatish Balay         sigma = 0.0;
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay         if (needH) {
228ffad9901SBarry Smith           ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
229a7e14dcfSSatish Balay           needH = 0;
230a7e14dcfSSatish Balay         }
231a7e14dcfSSatish Balay 
232a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
233a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
234a7e14dcfSSatish Balay           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
235a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
236a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
237a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
23887f595a5SBarry Smith           } else {
239a7e14dcfSSatish Balay             if (ftrial < fmin) {
240a7e14dcfSSatish Balay               fmin = ftrial;
241a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
242a7e14dcfSSatish Balay             }
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
245a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
248a7e14dcfSSatish Balay             actred = f - ftrial;
24987f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
250a7e14dcfSSatish Balay               kappa = 1.0;
25187f595a5SBarry Smith             } else {
252a7e14dcfSSatish Balay               kappa = actred / prered;
253a7e14dcfSSatish Balay             }
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
256a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
257a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
258a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
261a7e14dcfSSatish Balay               /* Great agreement */
262a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
263a7e14dcfSSatish Balay 
264a7e14dcfSSatish Balay               if (tau_max < 1.0) {
265a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
26687f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
267a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
26887f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
269a7e14dcfSSatish Balay                 tau = tau_1;
27087f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
271a7e14dcfSSatish Balay                 tau = tau_2;
27287f595a5SBarry Smith               } else {
273a7e14dcfSSatish Balay                 tau = tau_max;
274a7e14dcfSSatish Balay               }
27587f595a5SBarry Smith             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
276a7e14dcfSSatish Balay               /* Good agreement */
277a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
278a7e14dcfSSatish Balay 
279a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
280a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
28187f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
282a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
28387f595a5SBarry Smith               } else {
284a7e14dcfSSatish Balay                 tau = tau_max;
285a7e14dcfSSatish Balay               }
28687f595a5SBarry Smith             } else {
287a7e14dcfSSatish Balay               /* Not good agreement */
288a7e14dcfSSatish Balay               if (tau_min > 1.0) {
289a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
29087f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
291a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
29287f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
293a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
29487f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
295a7e14dcfSSatish Balay                 tau = tau_1;
29687f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
297a7e14dcfSSatish Balay                 tau = tau_2;
29887f595a5SBarry Smith               } else {
299a7e14dcfSSatish Balay                 tau = tau_max;
300a7e14dcfSSatish Balay               }
301a7e14dcfSSatish Balay             }
302a7e14dcfSSatish Balay           }
303a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
304a7e14dcfSSatish Balay         }
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay         if (fmin < f) {
307a7e14dcfSSatish Balay           f = fmin;
308a7e14dcfSSatish Balay           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
309a7e14dcfSSatish Balay           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
310a7e14dcfSSatish Balay 
311*a9603a14SPatrick Farrell           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
31287f595a5SBarry Smith           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
313a7e14dcfSSatish Balay           needH = 1;
314a7e14dcfSSatish Balay 
3158931d482SJason Sarich           ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
31687f595a5SBarry Smith           if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
317a7e14dcfSSatish Balay         }
318a7e14dcfSSatish Balay       }
319a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
320a7e14dcfSSatish Balay 
321a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
322a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
323a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
324a7e14dcfSSatish Balay       break;
325a7e14dcfSSatish Balay 
326a7e14dcfSSatish Balay     default:
327a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
328a7e14dcfSSatish Balay       tao->trust = 0.0;
329a7e14dcfSSatish Balay       break;
330a7e14dcfSSatish Balay     }
331a7e14dcfSSatish Balay   }
332a7e14dcfSSatish Balay 
333a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
334a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
335a7e14dcfSSatish Balay      since the function value may have decreased */
336a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type) {
337a7e14dcfSSatish Balay     if (f != 0.0) {
338a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
33987f595a5SBarry Smith     } else {
340a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
341a7e14dcfSSatish Balay     }
342a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
343a7e14dcfSSatish Balay   }
344a7e14dcfSSatish Balay 
345a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
346a7e14dcfSSatish Balay   nlsP->newt = 0;
347a7e14dcfSSatish Balay   nlsP->bfgs = 0;
348a7e14dcfSSatish Balay   nlsP->sgrad = 0;
349a7e14dcfSSatish Balay   nlsP->grad = 0;
350a7e14dcfSSatish Balay 
351a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
352a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
3538931d482SJason Sarich     ++tao->niter;
354ae93cb3cSJason Sarich     tao->ksp_its=0;
355a7e14dcfSSatish Balay 
356a7e14dcfSSatish Balay     /* Compute the Hessian */
357a7e14dcfSSatish Balay     if (needH) {
358ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
359a7e14dcfSSatish Balay       needH = 0;
360a7e14dcfSSatish Balay     }
361a7e14dcfSSatish Balay 
36287f595a5SBarry Smith     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
363a7e14dcfSSatish Balay       /* Obtain diagonal for the bfgs preconditioner  */
364a7e14dcfSSatish Balay       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
365a7e14dcfSSatish Balay       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
366a7e14dcfSSatish Balay       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
367a7e14dcfSSatish Balay       ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
368a7e14dcfSSatish Balay     }
369a7e14dcfSSatish Balay 
370a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
371a7e14dcfSSatish Balay     if (pert > 0) {
372302440fdSBarry Smith       ierr = MatShift(tao->hessian, pert);CHKERRQ(ierr);
373a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
374a7e14dcfSSatish Balay         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
375a7e14dcfSSatish Balay       }
376a7e14dcfSSatish Balay     }
377a7e14dcfSSatish Balay 
378a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type) {
379a7e14dcfSSatish Balay       if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) {
380a7e14dcfSSatish Balay         /* Obtain diagonal for the bfgs preconditioner  */
381a7e14dcfSSatish Balay         ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
382a7e14dcfSSatish Balay         ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
383a7e14dcfSSatish Balay         ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
384a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
385a7e14dcfSSatish Balay       }
386a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
387a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
388a7e14dcfSSatish Balay       ++bfgsUpdates;
389a7e14dcfSSatish Balay     }
390a7e14dcfSSatish Balay 
391a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
39223ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
39387f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type ||  NLS_KSP_GLTR == nlsP->ksp_type) {
394a7e14dcfSSatish Balay 
395a7e14dcfSSatish Balay       if (NLS_KSP_NASH == nlsP->ksp_type) {
396a7e14dcfSSatish Balay         ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
397a7e14dcfSSatish Balay       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
398a7e14dcfSSatish Balay          ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
399a7e14dcfSSatish Balay       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
400a7e14dcfSSatish Balay         ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
401a7e14dcfSSatish Balay       }
402a7e14dcfSSatish Balay 
403a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
404a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
405a7e14dcfSSatish Balay       tao->ksp_its+=kspits;
406ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
407a7e14dcfSSatish Balay 
408a7e14dcfSSatish Balay       if (NLS_KSP_NASH == nlsP->ksp_type) {
409a7e14dcfSSatish Balay         ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
410a7e14dcfSSatish Balay       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
411a7e14dcfSSatish Balay          ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
412a7e14dcfSSatish Balay       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
413a7e14dcfSSatish Balay         ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
414a7e14dcfSSatish Balay       }
415a7e14dcfSSatish Balay 
416a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
417a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
418a7e14dcfSSatish Balay         if (norm_d > 0.0) {
419a7e14dcfSSatish Balay           tao->trust = norm_d;
420a7e14dcfSSatish Balay 
421a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
422a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
423a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
42487f595a5SBarry Smith         } else {
425a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
426a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
427a7e14dcfSSatish Balay           tao->trust = tao->trust0;
428a7e14dcfSSatish Balay 
429a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
430a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
431a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
432a7e14dcfSSatish Balay 
433a7e14dcfSSatish Balay           if (NLS_KSP_NASH == nlsP->ksp_type) {
434a7e14dcfSSatish Balay             ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
435a7e14dcfSSatish Balay           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
436a7e14dcfSSatish Balay             ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
437a7e14dcfSSatish Balay           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
438a7e14dcfSSatish Balay             ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
439a7e14dcfSSatish Balay           }
440a7e14dcfSSatish Balay 
441a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
442a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
443a7e14dcfSSatish Balay           tao->ksp_its+=kspits;
444ae93cb3cSJason Sarich           tao->ksp_tot_its+=kspits;
445a7e14dcfSSatish Balay           if (NLS_KSP_NASH == nlsP->ksp_type) {
446a7e14dcfSSatish Balay             ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
447a7e14dcfSSatish Balay           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
448a7e14dcfSSatish Balay             ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
449a7e14dcfSSatish Balay           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
450a7e14dcfSSatish Balay             ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
451a7e14dcfSSatish Balay           }
452a7e14dcfSSatish Balay 
45387f595a5SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
454a7e14dcfSSatish Balay         }
455a7e14dcfSSatish Balay       }
45687f595a5SBarry Smith     } else {
457a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
458a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
459a7e14dcfSSatish Balay       tao->ksp_its += kspits;
460ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
461a7e14dcfSSatish Balay     }
462a7e14dcfSSatish Balay     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
463a7e14dcfSSatish Balay     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
46487f595a5SBarry Smith     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
465a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
466a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
467a7e14dcfSSatish Balay 
468a7e14dcfSSatish Balay       if (f != 0.0) {
469a7e14dcfSSatish Balay         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
47087f595a5SBarry Smith       } else {
471a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
472a7e14dcfSSatish Balay       }
473a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
474a7e14dcfSSatish Balay       ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
475a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
476a7e14dcfSSatish Balay       bfgsUpdates = 1;
477a7e14dcfSSatish Balay     }
478a7e14dcfSSatish Balay 
479a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
480a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
48187f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
482a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
48387f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
484a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
48587f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
486a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
48787f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
488a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
48987f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
490a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
49187f595a5SBarry Smith     } else {
492a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
493a7e14dcfSSatish Balay     }
494a7e14dcfSSatish Balay 
495a7e14dcfSSatish Balay     /* Check for success (descent direction) */
496a7e14dcfSSatish Balay     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
497a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
498a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
499a7e14dcfSSatish Balay          Update the perturbation for next time */
500a7e14dcfSSatish Balay       if (pert <= 0.0) {
501a7e14dcfSSatish Balay         /* Initialize the perturbation */
502a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
503a7e14dcfSSatish Balay         if (NLS_KSP_GLTR == nlsP->ksp_type) {
504a7e14dcfSSatish Balay           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
505a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
506a7e14dcfSSatish Balay         }
50787f595a5SBarry Smith       } else {
508a7e14dcfSSatish Balay         /* Increase the perturbation */
509a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
510a7e14dcfSSatish Balay       }
511a7e14dcfSSatish Balay 
512a7e14dcfSSatish Balay       if (NLS_PC_BFGS != nlsP->pc_type) {
513a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
514a7e14dcfSSatish Balay            Must use gradient direction in this case */
515a7e14dcfSSatish Balay         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
516a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
517a7e14dcfSSatish Balay         ++nlsP->grad;
518a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
51987f595a5SBarry Smith       } else {
520a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
521a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
522a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
523a7e14dcfSSatish Balay 
524a7e14dcfSSatish Balay         /* Check for success (descent direction) */
525a7e14dcfSSatish Balay         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
526a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
527a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
528a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
529a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
530a7e14dcfSSatish Balay              which is guaranteed to be descent */
531a7e14dcfSSatish Balay 
532a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
533a7e14dcfSSatish Balay 
534a7e14dcfSSatish Balay           if (f != 0.0) {
535a7e14dcfSSatish Balay             delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
53687f595a5SBarry Smith           } else {
537a7e14dcfSSatish Balay             delta = 2.0 / (gnorm*gnorm);
538a7e14dcfSSatish Balay           }
539a7e14dcfSSatish Balay           ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
540a7e14dcfSSatish Balay           ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
541a7e14dcfSSatish Balay           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
542a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
543a7e14dcfSSatish Balay           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
544a7e14dcfSSatish Balay 
545a7e14dcfSSatish Balay           bfgsUpdates = 1;
546a7e14dcfSSatish Balay           ++nlsP->sgrad;
547a7e14dcfSSatish Balay           stepType = NLS_SCALED_GRADIENT;
54887f595a5SBarry Smith         } else {
549a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
550a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
551a7e14dcfSSatish Balay             ++nlsP->sgrad;
552a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
55387f595a5SBarry Smith           } else {
554a7e14dcfSSatish Balay             ++nlsP->bfgs;
555a7e14dcfSSatish Balay             stepType = NLS_BFGS;
556a7e14dcfSSatish Balay           }
557a7e14dcfSSatish Balay         }
558a7e14dcfSSatish Balay       }
55987f595a5SBarry Smith     } else {
560a7e14dcfSSatish Balay       /* Computed Newton step is descent */
561a7e14dcfSSatish Balay       switch (ksp_reason) {
562a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
563a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
564a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
565a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
566a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
567a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
568a7e14dcfSSatish Balay         if (pert <= 0.0) {
569a7e14dcfSSatish Balay           /* Initialize the perturbation */
570a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
571a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
572a7e14dcfSSatish Balay             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
573a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
574a7e14dcfSSatish Balay           }
57587f595a5SBarry Smith         } else {
576a7e14dcfSSatish Balay           /* Increase the perturbation */
577a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
578a7e14dcfSSatish Balay         }
579a7e14dcfSSatish Balay         break;
580a7e14dcfSSatish Balay 
581a7e14dcfSSatish Balay       default:
582a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
583a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
584a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
585a7e14dcfSSatish Balay           pert = 0.0;
586a7e14dcfSSatish Balay         }
587a7e14dcfSSatish Balay         break;
588a7e14dcfSSatish Balay       }
589a7e14dcfSSatish Balay 
590a7e14dcfSSatish Balay       ++nlsP->newt;
591a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
592a7e14dcfSSatish Balay     }
593a7e14dcfSSatish Balay 
594a7e14dcfSSatish Balay     /* Perform the linesearch */
595a7e14dcfSSatish Balay     fold = f;
596a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
597a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
598a7e14dcfSSatish Balay 
599a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
600a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
601a7e14dcfSSatish Balay 
60287f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
603a7e14dcfSSatish Balay       f = fold;
604a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
605a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
606a7e14dcfSSatish Balay 
607a7e14dcfSSatish Balay       switch(stepType) {
608a7e14dcfSSatish Balay       case NLS_NEWTON:
609a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
610a7e14dcfSSatish Balay            Update the perturbation for next time */
611a7e14dcfSSatish Balay         if (pert <= 0.0) {
612a7e14dcfSSatish Balay           /* Initialize the perturbation */
613a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
614a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
615a7e14dcfSSatish Balay             ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
616a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
617a7e14dcfSSatish Balay           }
61887f595a5SBarry Smith         } else {
619a7e14dcfSSatish Balay           /* Increase the perturbation */
620a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
621a7e14dcfSSatish Balay         }
622a7e14dcfSSatish Balay 
623a7e14dcfSSatish Balay         if (NLS_PC_BFGS != nlsP->pc_type) {
624a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
625a7e14dcfSSatish Balay              Must use gradient direction in this case */
626a7e14dcfSSatish Balay           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
627a7e14dcfSSatish Balay           ++nlsP->grad;
628a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
62987f595a5SBarry Smith         } else {
630a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
631a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
632a7e14dcfSSatish Balay           /* Check for success (descent direction) */
633a7e14dcfSSatish Balay           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
634a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
635a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
636a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
637a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
638a7e14dcfSSatish Balay 
639a7e14dcfSSatish Balay             if (f != 0.0) {
640a7e14dcfSSatish Balay               delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
64187f595a5SBarry Smith             } else {
642a7e14dcfSSatish Balay               delta = 2.0 / (gnorm*gnorm);
643a7e14dcfSSatish Balay             }
644a7e14dcfSSatish Balay             ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
645a7e14dcfSSatish Balay             ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
646a7e14dcfSSatish Balay             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
647a7e14dcfSSatish Balay             ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
648a7e14dcfSSatish Balay 
649a7e14dcfSSatish Balay             bfgsUpdates = 1;
650a7e14dcfSSatish Balay             ++nlsP->sgrad;
651a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
65287f595a5SBarry Smith           } else {
653a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
654a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
655a7e14dcfSSatish Balay               ++nlsP->sgrad;
656a7e14dcfSSatish Balay               stepType = NLS_SCALED_GRADIENT;
65787f595a5SBarry Smith             } else {
658a7e14dcfSSatish Balay               ++nlsP->bfgs;
659a7e14dcfSSatish Balay               stepType = NLS_BFGS;
660a7e14dcfSSatish Balay             }
661a7e14dcfSSatish Balay           }
662a7e14dcfSSatish Balay         }
663a7e14dcfSSatish Balay         break;
664a7e14dcfSSatish Balay 
665a7e14dcfSSatish Balay       case NLS_BFGS:
666a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
667a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
668a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
669a7e14dcfSSatish Balay 
670a7e14dcfSSatish Balay         if (f != 0.0) {
671a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
67287f595a5SBarry Smith         } else {
673a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
674a7e14dcfSSatish Balay         }
675a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
676a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
677a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
678a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
679a7e14dcfSSatish Balay 
680a7e14dcfSSatish Balay         bfgsUpdates = 1;
681a7e14dcfSSatish Balay         ++nlsP->sgrad;
682a7e14dcfSSatish Balay         stepType = NLS_SCALED_GRADIENT;
683a7e14dcfSSatish Balay         break;
684a7e14dcfSSatish Balay 
685a7e14dcfSSatish Balay       case NLS_SCALED_GRADIENT:
686a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
687a7e14dcfSSatish Balay            The scaled gradient step did not produce a new iterate;
688a7e14dcfSSatish Balay            attemp to use the gradient direction.
689a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
690a7e14dcfSSatish Balay 
691a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr);
692a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr);
693a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
694a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
695a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
696a7e14dcfSSatish Balay 
697a7e14dcfSSatish Balay         bfgsUpdates = 1;
698a7e14dcfSSatish Balay         ++nlsP->grad;
699a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
700a7e14dcfSSatish Balay         break;
701a7e14dcfSSatish Balay       }
702a7e14dcfSSatish Balay       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
703a7e14dcfSSatish Balay 
704a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
705a7e14dcfSSatish Balay       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
706a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
707a7e14dcfSSatish Balay     }
708a7e14dcfSSatish Balay 
70987f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
710a7e14dcfSSatish Balay       /* Failed to find an improving point */
711a7e14dcfSSatish Balay       f = fold;
712a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
713a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
714a7e14dcfSSatish Balay       step = 0.0;
715a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
716a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
717a7e14dcfSSatish Balay       break;
718a7e14dcfSSatish Balay     }
719a7e14dcfSSatish Balay 
720a7e14dcfSSatish Balay     /* Update trust region radius */
72187f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
722a7e14dcfSSatish Balay       switch(nlsP->update_type) {
723a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
724a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
725a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
726a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
727a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
72887f595a5SBarry Smith           } else if (step < nlsP->nu2) {
729a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
730a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
73187f595a5SBarry Smith           } else if (step < nlsP->nu3) {
732a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
733a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
734a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
73587f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
736a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
737a7e14dcfSSatish Balay             }
73887f595a5SBarry Smith           } else if (step < nlsP->nu4) {
739a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
740a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
74187f595a5SBarry Smith           } else {
742a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
743a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
744a7e14dcfSSatish Balay           }
74587f595a5SBarry Smith         } else {
746a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
747a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
748a7e14dcfSSatish Balay         }
749a7e14dcfSSatish Balay         break;
750a7e14dcfSSatish Balay 
751a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
752a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
753a7e14dcfSSatish Balay           /*  Get predicted reduction */
754a7e14dcfSSatish Balay 
755a7e14dcfSSatish Balay           if (NLS_KSP_STCG == nlsP->ksp_type) {
756302440fdSBarry Smith             ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
757a7e14dcfSSatish Balay           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
758302440fdSBarry Smith             ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
759a7e14dcfSSatish Balay           } else {
760a7e14dcfSSatish Balay             ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
761a7e14dcfSSatish Balay           }
762a7e14dcfSSatish Balay 
763a7e14dcfSSatish Balay           if (prered >= 0.0) {
764a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
765a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
766a7e14dcfSSatish Balay             /*  be rejected! */
767a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
76887f595a5SBarry Smith           } else {
769a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
770a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
77187f595a5SBarry Smith             } else {
772a7e14dcfSSatish Balay               /*  Compute and actual reduction */
773a7e14dcfSSatish Balay               actred = fold - f_full;
774a7e14dcfSSatish Balay               prered = -prered;
775a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
776a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
777a7e14dcfSSatish Balay                 kappa = 1.0;
77887f595a5SBarry Smith               } else {
779a7e14dcfSSatish Balay                 kappa = actred / prered;
780a7e14dcfSSatish Balay               }
781a7e14dcfSSatish Balay 
782a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
783a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
784a7e14dcfSSatish Balay                 /*  Very bad step */
785a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
78687f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
787a7e14dcfSSatish Balay                 /*  Marginal bad step */
788a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
78987f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
790a7e14dcfSSatish Balay                 /*  Reasonable step */
791a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
792a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
79387f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
794a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
795a7e14dcfSSatish Balay                 }
79687f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
797a7e14dcfSSatish Balay                 /*  Good step */
798a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
79987f595a5SBarry Smith               } else {
800a7e14dcfSSatish Balay                 /*  Very good step */
801a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
802a7e14dcfSSatish Balay               }
803a7e14dcfSSatish Balay             }
804a7e14dcfSSatish Balay           }
80587f595a5SBarry Smith         } else {
806a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
807a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
808a7e14dcfSSatish Balay         }
809a7e14dcfSSatish Balay         break;
810a7e14dcfSSatish Balay 
811a7e14dcfSSatish Balay       default:
812a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
813a7e14dcfSSatish Balay 
814a7e14dcfSSatish Balay           if (NLS_KSP_STCG == nlsP->ksp_type) {
815302440fdSBarry Smith               ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
816a7e14dcfSSatish Balay           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
817302440fdSBarry Smith               ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
818a7e14dcfSSatish Balay           } else {
819a7e14dcfSSatish Balay               ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
820a7e14dcfSSatish Balay           }
821a7e14dcfSSatish Balay           if (prered >= 0.0) {
822a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
823a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
824a7e14dcfSSatish Balay             /*  be rejected! */
825a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
82687f595a5SBarry Smith           } else {
827a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
828a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
82987f595a5SBarry Smith             } else {
830a7e14dcfSSatish Balay               actred = fold - f_full;
831a7e14dcfSSatish Balay               prered = -prered;
83287f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
833a7e14dcfSSatish Balay                 kappa = 1.0;
83487f595a5SBarry Smith               } else {
835a7e14dcfSSatish Balay                 kappa = actred / prered;
836a7e14dcfSSatish Balay               }
837a7e14dcfSSatish Balay 
838a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
839a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
840a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
841a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
842a7e14dcfSSatish Balay 
843a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
844a7e14dcfSSatish Balay                 /*  Great agreement */
845a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
846a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
84787f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
848a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
84987f595a5SBarry Smith                 } else {
850a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
851a7e14dcfSSatish Balay                 }
85287f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
853a7e14dcfSSatish Balay                 /*  Good agreement */
854a7e14dcfSSatish Balay 
855a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
856a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
85787f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
858a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
85987f595a5SBarry Smith                 } else if (tau_max < 1.0) {
860a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
86187f595a5SBarry Smith                 } else {
862a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
863a7e14dcfSSatish Balay                 }
86487f595a5SBarry Smith               } else {
865a7e14dcfSSatish Balay                 /*  Not good agreement */
866a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
867a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
86887f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
869a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
87087f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
871a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
87287f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
873a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
87487f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
875a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
87687f595a5SBarry Smith                 } else {
877a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
878a7e14dcfSSatish Balay                 }
879a7e14dcfSSatish Balay               }
880a7e14dcfSSatish Balay             }
881a7e14dcfSSatish Balay           }
88287f595a5SBarry Smith         } else {
883a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
884a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
885a7e14dcfSSatish Balay         }
886a7e14dcfSSatish Balay         break;
887a7e14dcfSSatish Balay       }
888a7e14dcfSSatish Balay 
889a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
890a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
891a7e14dcfSSatish Balay     }
892a7e14dcfSSatish Balay 
893a7e14dcfSSatish Balay     /*  Check for termination */
894*a9603a14SPatrick Farrell     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
89587f595a5SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
896a7e14dcfSSatish Balay     needH = 1;
8978931d482SJason Sarich     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
898a7e14dcfSSatish Balay   }
899a7e14dcfSSatish Balay   PetscFunctionReturn(0);
900a7e14dcfSSatish Balay }
901a7e14dcfSSatish Balay 
902a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
903a7e14dcfSSatish Balay #undef __FUNCT__
904a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NLS"
905441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao)
906a7e14dcfSSatish Balay {
907a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
908a7e14dcfSSatish Balay   PetscErrorCode ierr;
909a7e14dcfSSatish Balay 
910a7e14dcfSSatish Balay   PetscFunctionBegin;
911a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
912a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
913a7e14dcfSSatish Balay   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
914a7e14dcfSSatish Balay   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
915a7e14dcfSSatish Balay   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
916a7e14dcfSSatish Balay   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
917a7e14dcfSSatish Balay   nlsP->Diag = 0;
918a7e14dcfSSatish Balay   nlsP->M = 0;
919a7e14dcfSSatish Balay   PetscFunctionReturn(0);
920a7e14dcfSSatish Balay }
921a7e14dcfSSatish Balay 
922a7e14dcfSSatish Balay /*------------------------------------------------------------*/
923a7e14dcfSSatish Balay #undef __FUNCT__
924a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NLS"
925441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao)
926a7e14dcfSSatish Balay {
927a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
928a7e14dcfSSatish Balay   PetscErrorCode ierr;
929a7e14dcfSSatish Balay 
930a7e14dcfSSatish Balay   PetscFunctionBegin;
931a7e14dcfSSatish Balay   if (tao->setupcalled) {
932a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
933a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
934a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
935a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
936a7e14dcfSSatish Balay   }
937a7e14dcfSSatish Balay   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
938a7e14dcfSSatish Balay   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
939a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
940a7e14dcfSSatish Balay   PetscFunctionReturn(0);
941a7e14dcfSSatish Balay }
942a7e14dcfSSatish Balay 
943a7e14dcfSSatish Balay /*------------------------------------------------------------*/
944a7e14dcfSSatish Balay #undef __FUNCT__
945a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NLS"
9468c34d3f5SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptions *PetscOptionsObject,Tao tao)
947a7e14dcfSSatish Balay {
948a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
949a7e14dcfSSatish Balay   PetscErrorCode ierr;
950a7e14dcfSSatish Balay 
951a7e14dcfSSatish Balay   PetscFunctionBegin;
9521a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
953a7e14dcfSSatish 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);
954a7e14dcfSSatish 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);
955a7e14dcfSSatish 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);
956a7e14dcfSSatish 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);
957a7e14dcfSSatish 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);
95894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr);
95994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr);
96094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr);
96194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr);
96294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr);
96394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr);
96494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr);
96594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr);
96694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr);
96794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr);
96894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr);
96994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr);
97094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr);
97194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr);
97294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr);
97394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr);
97494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr);
97594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr);
97694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr);
97794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr);
97894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr);
97994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr);
98094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr);
98194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr);
98294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr);
98394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr);
98494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr);
98594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr);
98694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr);
98794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr);
98894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr);
98994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr);
99094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr);
99194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr);
99294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr);
99394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr);
99494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr);
99594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr);
99694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr);
99794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr);
99894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr);
99994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr);
100094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr);
100194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr);
100294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr);
1003a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
1004a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1005a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1006a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1007a7e14dcfSSatish Balay }
1008a7e14dcfSSatish Balay 
1009a7e14dcfSSatish Balay 
1010a7e14dcfSSatish Balay /*------------------------------------------------------------*/
1011a7e14dcfSSatish Balay #undef __FUNCT__
1012a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NLS"
1013441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
1014a7e14dcfSSatish Balay {
1015a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
1016a7e14dcfSSatish Balay   PetscInt       nrejects;
1017a7e14dcfSSatish Balay   PetscBool      isascii;
1018a7e14dcfSSatish Balay   PetscErrorCode ierr;
1019a7e14dcfSSatish Balay 
1020a7e14dcfSSatish Balay   PetscFunctionBegin;
1021a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1022a7e14dcfSSatish Balay   if (isascii) {
1023a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1024a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
1025a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr);
1026a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1027a7e14dcfSSatish Balay     }
1028a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
1029a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
1030a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
1031a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
1032a7e14dcfSSatish Balay 
1033a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
1034a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
1035a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
1036a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
1037a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
1038a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
1039a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
1040a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1041a7e14dcfSSatish Balay   }
1042a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1043a7e14dcfSSatish Balay }
1044a7e14dcfSSatish Balay 
1045a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
10464aa34175SJason Sarich /*MC
10474aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
10484aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
10494aa34175SJason Sarich   system of equations to obtain the step diretion dk:
10504aa34175SJason Sarich               Hk dk = -gk
10514aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
10524aa34175SJason Sarich   solve
10534aa34175SJason Sarich               min_t f(xk + t d_k)
10544aa34175SJason Sarich 
10554aa34175SJason Sarich     Options Database Keys:
10564aa34175SJason Sarich + -tao_nls_ksp_type - "cg","nash","stcg","gltr","petsc"
10574aa34175SJason Sarich . -tao_nls_pc_type - "none","ahess","bfgs","petsc"
10584aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs"
10594aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation"
10604aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
10614aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
10624aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
10634aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
10644aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
10654aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
10664aa34175SJason Sarich . -tao_nls_pgfac - growth factor
10674aa34175SJason Sarich . -tao_nls_psfac - shrink factor
10684aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
10694aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
10704aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
10714aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
10724aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
10734aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
10744aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
10754aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
10764aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
10774aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
10784aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
10794aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
10804aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
10814aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
10824aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
10834aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
10844aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
10854aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
10864aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
10874aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
10884aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
10894aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
10904aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
10914aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
10921522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
10931522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
10941522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
10951522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
10961522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
10971522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
10981522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
10991eb8069cSJason Sarich 
11001eb8069cSJason Sarich   Level: beginner
11011522df2eSJason Sarich M*/
11024aa34175SJason Sarich 
1103a7e14dcfSSatish Balay #undef __FUNCT__
1104a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NLS"
1105728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
1106a7e14dcfSSatish Balay {
1107a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
11088caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
1109a7e14dcfSSatish Balay   PetscErrorCode ierr;
1110a7e14dcfSSatish Balay 
1111a7e14dcfSSatish Balay   PetscFunctionBegin;
11123c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1113a7e14dcfSSatish Balay 
1114a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
1115a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
1116a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
1117a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1118a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
1119a7e14dcfSSatish Balay 
11206552cf8aSJason Sarich   /* Override default settings (unless already changed) */
11216552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
11226552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
11236f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
11246552cf8aSJason Sarich   if (!tao->fatol_changed) tao->fatol = 1.0e-5;
11256552cf8aSJason Sarich   if (!tao->frtol_changed) tao->frtol = 1.0e-5;
11266f4723b1SBarry Smith #else
11276552cf8aSJason Sarich   if (!tao->fatol_changed) tao->fatol = 1.0e-10;
11286552cf8aSJason Sarich   if (!tao->frtol_changed) tao->frtol = 1.0e-10;
11296f4723b1SBarry Smith #endif
11306552cf8aSJason Sarich 
1131a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
1132a7e14dcfSSatish Balay 
1133a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
1134a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
1135a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
1136a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
1137a7e14dcfSSatish Balay 
1138a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
1139a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
1140a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
1141a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
1142a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
1143a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
1144a7e14dcfSSatish Balay 
1145a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
1146a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
1147a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
1148a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
1149a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
1150a7e14dcfSSatish Balay 
1151a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
1152a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
1153a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
1154a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
1155a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
1156a7e14dcfSSatish Balay 
1157a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
1158a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
1159a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
1160a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
1161a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
1162a7e14dcfSSatish Balay 
1163a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
1164a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
1165a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
1166a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
1167a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
1168a7e14dcfSSatish Balay 
1169a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
1170a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
1171a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
1172a7e14dcfSSatish Balay 
1173a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
1174a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
1175a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
1176a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
1177a7e14dcfSSatish Balay 
1178a7e14dcfSSatish Balay   nlsP->theta = 0.05;
1179a7e14dcfSSatish Balay 
1180a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
1181a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
1182a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
1183a7e14dcfSSatish Balay 
1184a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
1185a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
1186a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
1187a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
1188a7e14dcfSSatish Balay 
1189a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
1190a7e14dcfSSatish Balay 
1191a7e14dcfSSatish Balay   /*  Remaining parameters */
1192a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
1193a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
1194a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
1195a7e14dcfSSatish Balay 
1196a7e14dcfSSatish Balay   nlsP->ksp_type        = NLS_KSP_STCG;
1197a7e14dcfSSatish Balay   nlsP->pc_type         = NLS_PC_BFGS;
1198a7e14dcfSSatish Balay   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1199a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1200a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
1201a7e14dcfSSatish Balay 
1202a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1203a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1204441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
12055d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1206a7e14dcfSSatish Balay 
1207a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
1208a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
12095d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
1210a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1211a7e14dcfSSatish Balay }
1212a7e14dcfSSatish Balay 
1213a7e14dcfSSatish Balay #undef __FUNCT__
1214a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
1215a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
1216a7e14dcfSSatish Balay {
1217a7e14dcfSSatish Balay   PetscErrorCode ierr;
1218a7e14dcfSSatish Balay   Mat            M;
121987f595a5SBarry Smith 
1220a7e14dcfSSatish Balay   PetscFunctionBegin;
1221a7e14dcfSSatish Balay   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1222a7e14dcfSSatish Balay   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
1223a7e14dcfSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
1224a7e14dcfSSatish Balay   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
1225a7e14dcfSSatish Balay   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
1226a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1227a7e14dcfSSatish Balay }
1228