xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision ba7fe8fb864bc2e28543a5593f37aff2eb73db37)
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>
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay #define NLS_KSP_CG      0
8a7e14dcfSSatish Balay #define NLS_KSP_NASH    1
9a7e14dcfSSatish Balay #define NLS_KSP_STCG    2
10a7e14dcfSSatish Balay #define NLS_KSP_GLTR    3
11a7e14dcfSSatish Balay #define NLS_KSP_PETSC   4
12a7e14dcfSSatish Balay #define NLS_KSP_TYPES   5
13a7e14dcfSSatish Balay 
14a7e14dcfSSatish Balay #define NLS_PC_NONE     0
15a7e14dcfSSatish Balay #define NLS_PC_AHESS    1
16a7e14dcfSSatish Balay #define NLS_PC_BFGS     2
17a7e14dcfSSatish Balay #define NLS_PC_PETSC    3
18a7e14dcfSSatish Balay #define NLS_PC_TYPES    4
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS        0
21a7e14dcfSSatish Balay #define BFGS_SCALE_PHESS        1
22a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS         2
23a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES        3
24a7e14dcfSSatish Balay 
25a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT         0
26a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION        1
27a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION    2
28a7e14dcfSSatish Balay #define NLS_INIT_TYPES            3
29a7e14dcfSSatish Balay 
30a7e14dcfSSatish Balay #define NLS_UPDATE_STEP           0
31a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION      1
32a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION  2
33a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES          3
34a7e14dcfSSatish Balay 
3587f595a5SBarry Smith static const char *NLS_KSP[64] = {"cg", "nash", "stcg", "gltr", "petsc"};
36a7e14dcfSSatish Balay 
3787f595a5SBarry Smith static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"};
38a7e14dcfSSatish Balay 
3987f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
40a7e14dcfSSatish Balay 
4187f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"};
42a7e14dcfSSatish Balay 
4387f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"};
44a7e14dcfSSatish Balay 
45a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x);
46a7e14dcfSSatish Balay /* Routine for BFGS preconditioner
47a7e14dcfSSatish Balay 
48a7e14dcfSSatish Balay 
49a7e14dcfSSatish Balay  Implements Newton's Method with a line search approach for solving
50a7e14dcfSSatish Balay  unconstrained minimization problems.  A More'-Thuente line search
51a7e14dcfSSatish Balay  is used to guarantee that the bfgs preconditioner remains positive
52a7e14dcfSSatish Balay  definite.
53a7e14dcfSSatish Balay 
54a7e14dcfSSatish Balay  The method can shift the Hessian matrix.  The shifting procedure is
55a7e14dcfSSatish Balay  adapted from the PATH algorithm for solving complementarity
56a7e14dcfSSatish Balay  problems.
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay  The linear system solve should be done with a conjugate gradient
59a7e14dcfSSatish Balay  method, although any method can be used. */
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay #define NLS_NEWTON              0
62a7e14dcfSSatish Balay #define NLS_BFGS                1
63a7e14dcfSSatish Balay #define NLS_SCALED_GRADIENT     2
64a7e14dcfSSatish Balay #define NLS_GRADIENT            3
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay #undef __FUNCT__
67a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NLS"
68441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao tao)
69a7e14dcfSSatish Balay {
70a7e14dcfSSatish Balay   PetscErrorCode               ierr;
71a7e14dcfSSatish Balay   TAO_NLS                      *nlsP = (TAO_NLS *)tao->data;
72a7e14dcfSSatish Balay   PC                           pc;
73a7e14dcfSSatish Balay   KSPConvergedReason           ksp_reason;
74e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
75e4cb33bbSBarry Smith   TaoConvergedReason           reason;
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay   PetscReal                    fmin, ftrial, f_full, prered, actred, kappa, sigma;
78a7e14dcfSSatish Balay   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
79a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm, pert;
80a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
81a7e14dcfSSatish Balay   PetscReal                    delta;
82a7e14dcfSSatish Balay   PetscReal                    norm_d = 0.0, e_min;
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay   PetscInt                     stepType;
85a7e14dcfSSatish Balay   PetscInt                     bfgsUpdates = 0;
86a7e14dcfSSatish Balay   PetscInt                     n,N,kspits;
87b4690188SBarry Smith   PetscInt                     needH = 1;
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay   PetscInt                     i_max = 5;
90a7e14dcfSSatish Balay   PetscInt                     j_max = 1;
91a7e14dcfSSatish Balay   PetscInt                     i, j;
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay   PetscFunctionBegin;
94a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
95a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
96a7e14dcfSSatish Balay   }
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   /* Initialized variables */
99a7e14dcfSSatish Balay   pert = nlsP->sval;
100a7e14dcfSSatish Balay 
1012d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
102a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
103a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
104a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
105a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
106a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
107a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
108a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay   /* Modify the linear solver to a trust region method if desired */
111a7e14dcfSSatish Balay   switch(nlsP->ksp_type) {
112a7e14dcfSSatish Balay   case NLS_KSP_CG:
113a7e14dcfSSatish Balay     ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr);
11472b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
115a7e14dcfSSatish Balay     break;
116a7e14dcfSSatish Balay 
117a7e14dcfSSatish Balay   case NLS_KSP_NASH:
118*ba7fe8fbSTodd Munson     ierr = KSPSetType(tao->ksp, KSPCGNASH);CHKERRQ(ierr);
11972b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
120a7e14dcfSSatish Balay     break;
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay   case NLS_KSP_STCG:
123*ba7fe8fbSTodd Munson     ierr = KSPSetType(tao->ksp, KSPCGSTCG);CHKERRQ(ierr);
12472b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
125a7e14dcfSSatish Balay     break;
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay   case NLS_KSP_GLTR:
128*ba7fe8fbSTodd Munson     ierr = KSPSetType(tao->ksp, KSPCGGLTR);CHKERRQ(ierr);
12972b7fd4bSBarry Smith     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
130a7e14dcfSSatish Balay     break;
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay   default:
133a7e14dcfSSatish Balay     /* Use the method set by the ksp_type */
134a7e14dcfSSatish Balay     break;
135a7e14dcfSSatish Balay   }
136a7e14dcfSSatish Balay 
137a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
138*ba7fe8fbSTodd Munson      Command automatically ignored for other methods
139*ba7fe8fbSTodd Munson      Will be reset during the first iteration
140*ba7fe8fbSTodd Munson   */
141*ba7fe8fbSTodd Munson   ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
142a7e14dcfSSatish Balay 
14387f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
144a7e14dcfSSatish Balay     tao->trust = tao->trust0;
14587f595a5SBarry Smith     if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative");
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay     /* Modify the radius if it is too large or small */
148a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
149a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
150a7e14dcfSSatish Balay   }
151a7e14dcfSSatish Balay 
152a7e14dcfSSatish Balay   /* Get vectors we will need */
153a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
154a7e14dcfSSatish Balay     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
155a7e14dcfSSatish Balay     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
156a7e14dcfSSatish Balay     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
157a7e14dcfSSatish Balay     ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr);
158a7e14dcfSSatish Balay   }
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay   /* Check convergence criteria */
161a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
162a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
16387f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
164a7e14dcfSSatish Balay 
1658931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
16687f595a5SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   /* create vectors for the limited memory preconditioner */
16987f595a5SBarry Smith   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
170a7e14dcfSSatish Balay     if (!nlsP->Diag) {
171a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
172a7e14dcfSSatish Balay     }
173a7e14dcfSSatish Balay   }
174a7e14dcfSSatish Balay 
175a7e14dcfSSatish Balay   /* Modify the preconditioner to use the bfgs approximation */
176a7e14dcfSSatish Balay   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
177a7e14dcfSSatish Balay   switch(nlsP->pc_type) {
178a7e14dcfSSatish Balay   case NLS_PC_NONE:
179a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
1801a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
181a7e14dcfSSatish Balay     break;
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay   case NLS_PC_AHESS:
184a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
1851a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
186baa89ecbSBarry Smith     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
187a7e14dcfSSatish Balay     break;
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay   case NLS_PC_BFGS:
190a7e14dcfSSatish Balay     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
1911a1499c8SBarry Smith     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
192a7e14dcfSSatish Balay     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
193a7e14dcfSSatish Balay     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
194a7e14dcfSSatish Balay     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
195a7e14dcfSSatish Balay     break;
196a7e14dcfSSatish Balay 
197a7e14dcfSSatish Balay   default:
198a7e14dcfSSatish Balay     /* Use the pc method set by pc_type */
199a7e14dcfSSatish Balay     break;
200a7e14dcfSSatish Balay   }
201a7e14dcfSSatish Balay 
202a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
203a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
20487f595a5SBarry Smith   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
205a7e14dcfSSatish Balay     switch(nlsP->init_type) {
206a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
207a7e14dcfSSatish Balay       /* Use the initial radius specified */
208a7e14dcfSSatish Balay       break;
209a7e14dcfSSatish Balay 
210a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
211a7e14dcfSSatish Balay       /* Use the initial radius specified */
212a7e14dcfSSatish Balay       max_radius = 0.0;
213a7e14dcfSSatish Balay 
214a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
215a7e14dcfSSatish Balay         fmin = f;
216a7e14dcfSSatish Balay         sigma = 0.0;
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay         if (needH) {
219ffad9901SBarry Smith           ierr  = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
220a7e14dcfSSatish Balay           needH = 0;
221a7e14dcfSSatish Balay         }
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
224a7e14dcfSSatish Balay           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
225a7e14dcfSSatish Balay           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
226a7e14dcfSSatish Balay           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
227a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
228a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
22987f595a5SBarry Smith           } else {
230a7e14dcfSSatish Balay             if (ftrial < fmin) {
231a7e14dcfSSatish Balay               fmin = ftrial;
232a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
233a7e14dcfSSatish Balay             }
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
236a7e14dcfSSatish Balay             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
239a7e14dcfSSatish Balay             actred = f - ftrial;
24087f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
241a7e14dcfSSatish Balay               kappa = 1.0;
24287f595a5SBarry Smith             } else {
243a7e14dcfSSatish Balay               kappa = actred / prered;
244a7e14dcfSSatish Balay             }
245a7e14dcfSSatish Balay 
246a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
247a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
248a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
249a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
250a7e14dcfSSatish Balay 
251a7e14dcfSSatish Balay             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
252a7e14dcfSSatish Balay               /* Great agreement */
253a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay               if (tau_max < 1.0) {
256a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
25787f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
258a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
25987f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
260a7e14dcfSSatish Balay                 tau = tau_1;
26187f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
262a7e14dcfSSatish Balay                 tau = tau_2;
26387f595a5SBarry Smith               } else {
264a7e14dcfSSatish Balay                 tau = tau_max;
265a7e14dcfSSatish Balay               }
26687f595a5SBarry Smith             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
267a7e14dcfSSatish Balay               /* Good agreement */
268a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
269a7e14dcfSSatish Balay 
270a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
271a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
27287f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
273a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
27487f595a5SBarry Smith               } else {
275a7e14dcfSSatish Balay                 tau = tau_max;
276a7e14dcfSSatish Balay               }
27787f595a5SBarry Smith             } else {
278a7e14dcfSSatish Balay               /* Not good agreement */
279a7e14dcfSSatish Balay               if (tau_min > 1.0) {
280a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
28187f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
282a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
28387f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
284a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
28587f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
286a7e14dcfSSatish Balay                 tau = tau_1;
28787f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
288a7e14dcfSSatish Balay                 tau = tau_2;
28987f595a5SBarry Smith               } else {
290a7e14dcfSSatish Balay                 tau = tau_max;
291a7e14dcfSSatish Balay               }
292a7e14dcfSSatish Balay             }
293a7e14dcfSSatish Balay           }
294a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
295a7e14dcfSSatish Balay         }
296a7e14dcfSSatish Balay 
297a7e14dcfSSatish Balay         if (fmin < f) {
298a7e14dcfSSatish Balay           f = fmin;
299a7e14dcfSSatish Balay           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
300a7e14dcfSSatish Balay           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
301a7e14dcfSSatish Balay 
302a9603a14SPatrick Farrell           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
30387f595a5SBarry Smith           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
304a7e14dcfSSatish Balay           needH = 1;
305a7e14dcfSSatish Balay 
3068931d482SJason Sarich           ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
30787f595a5SBarry Smith           if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
308a7e14dcfSSatish Balay         }
309a7e14dcfSSatish Balay       }
310a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
311a7e14dcfSSatish Balay 
312a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
313a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
314a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
315a7e14dcfSSatish Balay       break;
316a7e14dcfSSatish Balay 
317a7e14dcfSSatish Balay     default:
318a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
319a7e14dcfSSatish Balay       tao->trust = 0.0;
320a7e14dcfSSatish Balay       break;
321a7e14dcfSSatish Balay     }
322a7e14dcfSSatish Balay   }
323a7e14dcfSSatish Balay 
324a7e14dcfSSatish Balay   /* Set initial scaling for the BFGS preconditioner
325a7e14dcfSSatish Balay      This step is done after computing the initial trust-region radius
326a7e14dcfSSatish Balay      since the function value may have decreased */
327a7e14dcfSSatish Balay   if (NLS_PC_BFGS == nlsP->pc_type) {
328a7e14dcfSSatish Balay     if (f != 0.0) {
329a7e14dcfSSatish Balay       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
33087f595a5SBarry Smith     } else {
331a7e14dcfSSatish Balay       delta = 2.0 / (gnorm*gnorm);
332a7e14dcfSSatish Balay     }
333a7e14dcfSSatish Balay     ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
334a7e14dcfSSatish Balay   }
335a7e14dcfSSatish Balay 
336a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
337a7e14dcfSSatish Balay   nlsP->newt = 0;
338a7e14dcfSSatish Balay   nlsP->bfgs = 0;
339a7e14dcfSSatish Balay   nlsP->sgrad = 0;
340a7e14dcfSSatish Balay   nlsP->grad = 0;
341a7e14dcfSSatish Balay 
342a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
343a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
3448931d482SJason Sarich     ++tao->niter;
345ae93cb3cSJason Sarich     tao->ksp_its=0;
346a7e14dcfSSatish Balay 
347a7e14dcfSSatish Balay     /* Compute the Hessian */
348a7e14dcfSSatish Balay     if (needH) {
349ffad9901SBarry Smith       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
350a7e14dcfSSatish Balay     }
351a7e14dcfSSatish Balay 
35287f595a5SBarry Smith     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
353a7e14dcfSSatish Balay       /* Obtain diagonal for the bfgs preconditioner  */
354a7e14dcfSSatish Balay       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
355a7e14dcfSSatish Balay       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
356a7e14dcfSSatish Balay       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
357a7e14dcfSSatish Balay       ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
358a7e14dcfSSatish Balay     }
359a7e14dcfSSatish Balay 
360a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
361a7e14dcfSSatish Balay     if (pert > 0) {
362302440fdSBarry Smith       ierr = MatShift(tao->hessian, pert);CHKERRQ(ierr);
363a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
364a7e14dcfSSatish Balay         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
365a7e14dcfSSatish Balay       }
366a7e14dcfSSatish Balay     }
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type) {
369a7e14dcfSSatish Balay       if (BFGS_SCALE_PHESS == 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       /* Update the limited memory preconditioner */
377a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
378a7e14dcfSSatish Balay       ++bfgsUpdates;
379a7e14dcfSSatish Balay     }
380a7e14dcfSSatish Balay 
381a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
38223ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
38387f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type ||  NLS_KSP_GLTR == nlsP->ksp_type) {
384*ba7fe8fbSTodd Munson       ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
385a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
386a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
387a7e14dcfSSatish Balay       tao->ksp_its+=kspits;
388ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
389*ba7fe8fbSTodd Munson       ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
390a7e14dcfSSatish Balay 
391a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
392a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
393a7e14dcfSSatish Balay         if (norm_d > 0.0) {
394a7e14dcfSSatish Balay           tao->trust = norm_d;
395a7e14dcfSSatish Balay 
396a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
397a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
398a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
39987f595a5SBarry Smith         } else {
400a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
401a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
402a7e14dcfSSatish Balay           tao->trust = tao->trust0;
403a7e14dcfSSatish Balay 
404a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
405a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
406a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
407a7e14dcfSSatish Balay 
408*ba7fe8fbSTodd Munson           ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
409a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
410a7e14dcfSSatish Balay           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
411a7e14dcfSSatish Balay           tao->ksp_its+=kspits;
412ae93cb3cSJason Sarich           tao->ksp_tot_its+=kspits;
413*ba7fe8fbSTodd Munson           ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
414a7e14dcfSSatish Balay 
41587f595a5SBarry Smith           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
416a7e14dcfSSatish Balay         }
417a7e14dcfSSatish Balay       }
41887f595a5SBarry Smith     } else {
419a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
420a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
421a7e14dcfSSatish Balay       tao->ksp_its += kspits;
422ae93cb3cSJason Sarich       tao->ksp_tot_its+=kspits;
423a7e14dcfSSatish Balay     }
424a7e14dcfSSatish Balay     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
425a7e14dcfSSatish Balay     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
42687f595a5SBarry Smith     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
427a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
428a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay       if (f != 0.0) {
431a7e14dcfSSatish Balay         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
43287f595a5SBarry Smith       } else {
433a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
434a7e14dcfSSatish Balay       }
435a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
436a7e14dcfSSatish Balay       ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
437a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
438a7e14dcfSSatish Balay       bfgsUpdates = 1;
439a7e14dcfSSatish Balay     }
440a7e14dcfSSatish Balay 
441a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
442a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
44387f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
444a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
44587f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
446a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
44787f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
448a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
44987f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
450a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
45187f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
452a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
45387f595a5SBarry Smith     } else {
454a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
455a7e14dcfSSatish Balay     }
456a7e14dcfSSatish Balay 
457a7e14dcfSSatish Balay     /* Check for success (descent direction) */
458a7e14dcfSSatish Balay     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
459a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
460a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
461a7e14dcfSSatish Balay          Update the perturbation for next time */
462a7e14dcfSSatish Balay       if (pert <= 0.0) {
463a7e14dcfSSatish Balay         /* Initialize the perturbation */
464a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
465a7e14dcfSSatish Balay         if (NLS_KSP_GLTR == nlsP->ksp_type) {
466*ba7fe8fbSTodd Munson           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
467a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
468a7e14dcfSSatish Balay         }
46987f595a5SBarry Smith       } else {
470a7e14dcfSSatish Balay         /* Increase the perturbation */
471a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
472a7e14dcfSSatish Balay       }
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay       if (NLS_PC_BFGS != nlsP->pc_type) {
475a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
476a7e14dcfSSatish Balay            Must use gradient direction in this case */
477a7e14dcfSSatish Balay         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
478a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
479a7e14dcfSSatish Balay         ++nlsP->grad;
480a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
48187f595a5SBarry Smith       } else {
482a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
483a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
484a7e14dcfSSatish Balay         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
485a7e14dcfSSatish Balay 
486a7e14dcfSSatish Balay         /* Check for success (descent direction) */
487a7e14dcfSSatish Balay         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
488a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
489a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
490a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
491a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
492a7e14dcfSSatish Balay              which is guaranteed to be descent */
493a7e14dcfSSatish Balay 
494a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
495a7e14dcfSSatish Balay 
496a7e14dcfSSatish Balay           if (f != 0.0) {
497a7e14dcfSSatish Balay             delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
49887f595a5SBarry Smith           } else {
499a7e14dcfSSatish Balay             delta = 2.0 / (gnorm*gnorm);
500a7e14dcfSSatish Balay           }
501a7e14dcfSSatish Balay           ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
502a7e14dcfSSatish Balay           ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
503a7e14dcfSSatish Balay           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
504a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
505a7e14dcfSSatish Balay           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
506a7e14dcfSSatish Balay 
507a7e14dcfSSatish Balay           bfgsUpdates = 1;
508a7e14dcfSSatish Balay           ++nlsP->sgrad;
509a7e14dcfSSatish Balay           stepType = NLS_SCALED_GRADIENT;
51087f595a5SBarry Smith         } else {
511a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
512a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
513a7e14dcfSSatish Balay             ++nlsP->sgrad;
514a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
51587f595a5SBarry Smith           } else {
516a7e14dcfSSatish Balay             ++nlsP->bfgs;
517a7e14dcfSSatish Balay             stepType = NLS_BFGS;
518a7e14dcfSSatish Balay           }
519a7e14dcfSSatish Balay         }
520a7e14dcfSSatish Balay       }
52187f595a5SBarry Smith     } else {
522a7e14dcfSSatish Balay       /* Computed Newton step is descent */
523a7e14dcfSSatish Balay       switch (ksp_reason) {
524a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
525a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
526a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
527a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
528a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
529a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
530a7e14dcfSSatish Balay         if (pert <= 0.0) {
531a7e14dcfSSatish Balay           /* Initialize the perturbation */
532a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
533a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
534*ba7fe8fbSTodd Munson             ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
535a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
536a7e14dcfSSatish Balay           }
53787f595a5SBarry Smith         } else {
538a7e14dcfSSatish Balay           /* Increase the perturbation */
539a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
540a7e14dcfSSatish Balay         }
541a7e14dcfSSatish Balay         break;
542a7e14dcfSSatish Balay 
543a7e14dcfSSatish Balay       default:
544a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
545a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
546a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
547a7e14dcfSSatish Balay           pert = 0.0;
548a7e14dcfSSatish Balay         }
549a7e14dcfSSatish Balay         break;
550a7e14dcfSSatish Balay       }
551a7e14dcfSSatish Balay 
552a7e14dcfSSatish Balay       ++nlsP->newt;
553a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
554a7e14dcfSSatish Balay     }
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay     /* Perform the linesearch */
557a7e14dcfSSatish Balay     fold = f;
558a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
559a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
560a7e14dcfSSatish Balay 
561a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
562a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
563a7e14dcfSSatish Balay 
56487f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
565a7e14dcfSSatish Balay       f = fold;
566a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
567a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
568a7e14dcfSSatish Balay 
569a7e14dcfSSatish Balay       switch(stepType) {
570a7e14dcfSSatish Balay       case NLS_NEWTON:
571a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
572a7e14dcfSSatish Balay            Update the perturbation for next time */
573a7e14dcfSSatish Balay         if (pert <= 0.0) {
574a7e14dcfSSatish Balay           /* Initialize the perturbation */
575a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
576a7e14dcfSSatish Balay           if (NLS_KSP_GLTR == nlsP->ksp_type) {
577*ba7fe8fbSTodd Munson             ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
578a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
579a7e14dcfSSatish Balay           }
58087f595a5SBarry Smith         } else {
581a7e14dcfSSatish Balay           /* Increase the perturbation */
582a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
583a7e14dcfSSatish Balay         }
584a7e14dcfSSatish Balay 
585a7e14dcfSSatish Balay         if (NLS_PC_BFGS != nlsP->pc_type) {
586a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
587a7e14dcfSSatish Balay              Must use gradient direction in this case */
588a7e14dcfSSatish Balay           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
589a7e14dcfSSatish Balay           ++nlsP->grad;
590a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
59187f595a5SBarry Smith         } else {
592a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
593a7e14dcfSSatish Balay           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
594a7e14dcfSSatish Balay           /* Check for success (descent direction) */
595a7e14dcfSSatish Balay           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
596a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
597a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
598a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
599a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay             if (f != 0.0) {
602a7e14dcfSSatish Balay               delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
60387f595a5SBarry Smith             } else {
604a7e14dcfSSatish Balay               delta = 2.0 / (gnorm*gnorm);
605a7e14dcfSSatish Balay             }
606a7e14dcfSSatish Balay             ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
607a7e14dcfSSatish Balay             ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
608a7e14dcfSSatish Balay             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
609a7e14dcfSSatish Balay             ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
610a7e14dcfSSatish Balay 
611a7e14dcfSSatish Balay             bfgsUpdates = 1;
612a7e14dcfSSatish Balay             ++nlsP->sgrad;
613a7e14dcfSSatish Balay             stepType = NLS_SCALED_GRADIENT;
61487f595a5SBarry Smith           } else {
615a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
616a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
617a7e14dcfSSatish Balay               ++nlsP->sgrad;
618a7e14dcfSSatish Balay               stepType = NLS_SCALED_GRADIENT;
61987f595a5SBarry Smith             } else {
620a7e14dcfSSatish Balay               ++nlsP->bfgs;
621a7e14dcfSSatish Balay               stepType = NLS_BFGS;
622a7e14dcfSSatish Balay             }
623a7e14dcfSSatish Balay           }
624a7e14dcfSSatish Balay         }
625a7e14dcfSSatish Balay         break;
626a7e14dcfSSatish Balay 
627a7e14dcfSSatish Balay       case NLS_BFGS:
628a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
629a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
630a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
631a7e14dcfSSatish Balay 
632a7e14dcfSSatish Balay         if (f != 0.0) {
633a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
63487f595a5SBarry Smith         } else {
635a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
636a7e14dcfSSatish Balay         }
637a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
638a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
639a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
640a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
641a7e14dcfSSatish Balay 
642a7e14dcfSSatish Balay         bfgsUpdates = 1;
643a7e14dcfSSatish Balay         ++nlsP->sgrad;
644a7e14dcfSSatish Balay         stepType = NLS_SCALED_GRADIENT;
645a7e14dcfSSatish Balay         break;
646a7e14dcfSSatish Balay 
647a7e14dcfSSatish Balay       case NLS_SCALED_GRADIENT:
648a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
649a7e14dcfSSatish Balay            The scaled gradient step did not produce a new iterate;
650a7e14dcfSSatish Balay            attemp to use the gradient direction.
651a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
652a7e14dcfSSatish Balay 
653a7e14dcfSSatish Balay         ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr);
654a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr);
655a7e14dcfSSatish Balay         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
656a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
657a7e14dcfSSatish Balay         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
658a7e14dcfSSatish Balay 
659a7e14dcfSSatish Balay         bfgsUpdates = 1;
660a7e14dcfSSatish Balay         ++nlsP->grad;
661a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
662a7e14dcfSSatish Balay         break;
663a7e14dcfSSatish Balay       }
664a7e14dcfSSatish Balay       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
665a7e14dcfSSatish Balay 
666a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
667a7e14dcfSSatish Balay       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
668a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
669a7e14dcfSSatish Balay     }
670a7e14dcfSSatish Balay 
67187f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
672a7e14dcfSSatish Balay       /* Failed to find an improving point */
673a7e14dcfSSatish Balay       f = fold;
674a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
675a7e14dcfSSatish Balay       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
676a7e14dcfSSatish Balay       step = 0.0;
677a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
678a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
679a7e14dcfSSatish Balay       break;
680a7e14dcfSSatish Balay     }
681a7e14dcfSSatish Balay 
682a7e14dcfSSatish Balay     /* Update trust region radius */
68387f595a5SBarry Smith     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
684a7e14dcfSSatish Balay       switch(nlsP->update_type) {
685a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
686a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
687a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
688a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
689a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
69087f595a5SBarry Smith           } else if (step < nlsP->nu2) {
691a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
692a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
69387f595a5SBarry Smith           } else if (step < nlsP->nu3) {
694a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
695a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
696a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
69787f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
698a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
699a7e14dcfSSatish Balay             }
70087f595a5SBarry Smith           } else if (step < nlsP->nu4) {
701a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
702a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
70387f595a5SBarry Smith           } else {
704a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
705a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
706a7e14dcfSSatish Balay           }
70787f595a5SBarry Smith         } else {
708a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
709a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
710a7e14dcfSSatish Balay         }
711a7e14dcfSSatish Balay         break;
712a7e14dcfSSatish Balay 
713a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
714a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
715a7e14dcfSSatish Balay           /*  Get predicted reduction */
716a7e14dcfSSatish Balay 
717*ba7fe8fbSTodd Munson           ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
718a7e14dcfSSatish Balay           if (prered >= 0.0) {
719a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
720a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
721a7e14dcfSSatish Balay             /*  be rejected! */
722a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
72387f595a5SBarry Smith           } else {
724a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
725a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
72687f595a5SBarry Smith             } else {
727a7e14dcfSSatish Balay               /*  Compute and actual reduction */
728a7e14dcfSSatish Balay               actred = fold - f_full;
729a7e14dcfSSatish Balay               prered = -prered;
730a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
731a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
732a7e14dcfSSatish Balay                 kappa = 1.0;
73387f595a5SBarry Smith               } else {
734a7e14dcfSSatish Balay                 kappa = actred / prered;
735a7e14dcfSSatish Balay               }
736a7e14dcfSSatish Balay 
737a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
738a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
739a7e14dcfSSatish Balay                 /*  Very bad step */
740a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
74187f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
742a7e14dcfSSatish Balay                 /*  Marginal bad step */
743a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
74487f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
745a7e14dcfSSatish Balay                 /*  Reasonable step */
746a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
747a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
74887f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
749a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
750a7e14dcfSSatish Balay                 }
75187f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
752a7e14dcfSSatish Balay                 /*  Good step */
753a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
75487f595a5SBarry Smith               } else {
755a7e14dcfSSatish Balay                 /*  Very good step */
756a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
757a7e14dcfSSatish Balay               }
758a7e14dcfSSatish Balay             }
759a7e14dcfSSatish Balay           }
76087f595a5SBarry Smith         } else {
761a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
762a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
763a7e14dcfSSatish Balay         }
764a7e14dcfSSatish Balay         break;
765a7e14dcfSSatish Balay 
766a7e14dcfSSatish Balay       default:
767a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
768*ba7fe8fbSTodd Munson           ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
769a7e14dcfSSatish Balay           if (prered >= 0.0) {
770a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
771a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
772a7e14dcfSSatish Balay             /*  be rejected! */
773a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
77487f595a5SBarry Smith           } else {
775a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
776a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
77787f595a5SBarry Smith             } else {
778a7e14dcfSSatish Balay               actred = fold - f_full;
779a7e14dcfSSatish Balay               prered = -prered;
78087f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
781a7e14dcfSSatish Balay                 kappa = 1.0;
78287f595a5SBarry Smith               } else {
783a7e14dcfSSatish Balay                 kappa = actred / prered;
784a7e14dcfSSatish Balay               }
785a7e14dcfSSatish Balay 
786a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
787a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
788a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
789a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
790a7e14dcfSSatish Balay 
791a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
792a7e14dcfSSatish Balay                 /*  Great agreement */
793a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
794a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
79587f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
796a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
79787f595a5SBarry Smith                 } else {
798a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
799a7e14dcfSSatish Balay                 }
80087f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
801a7e14dcfSSatish Balay                 /*  Good agreement */
802a7e14dcfSSatish Balay 
803a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
804a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
80587f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
806a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
80787f595a5SBarry Smith                 } else if (tau_max < 1.0) {
808a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
80987f595a5SBarry Smith                 } else {
810a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
811a7e14dcfSSatish Balay                 }
81287f595a5SBarry Smith               } else {
813a7e14dcfSSatish Balay                 /*  Not good agreement */
814a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
815a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
81687f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
817a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
81887f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
819a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
82087f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
821a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
82287f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
823a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
82487f595a5SBarry Smith                 } else {
825a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
826a7e14dcfSSatish Balay                 }
827a7e14dcfSSatish Balay               }
828a7e14dcfSSatish Balay             }
829a7e14dcfSSatish Balay           }
83087f595a5SBarry Smith         } else {
831a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
832a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
833a7e14dcfSSatish Balay         }
834a7e14dcfSSatish Balay         break;
835a7e14dcfSSatish Balay       }
836a7e14dcfSSatish Balay 
837a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
838a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
839a7e14dcfSSatish Balay     }
840a7e14dcfSSatish Balay 
841a7e14dcfSSatish Balay     /*  Check for termination */
842a9603a14SPatrick Farrell     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
84387f595a5SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
844a7e14dcfSSatish Balay     needH = 1;
8458931d482SJason Sarich     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
846a7e14dcfSSatish Balay   }
847a7e14dcfSSatish Balay   PetscFunctionReturn(0);
848a7e14dcfSSatish Balay }
849a7e14dcfSSatish Balay 
850a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
851a7e14dcfSSatish Balay #undef __FUNCT__
852a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NLS"
853441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao)
854a7e14dcfSSatish Balay {
855a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
856a7e14dcfSSatish Balay   PetscErrorCode ierr;
857a7e14dcfSSatish Balay 
858a7e14dcfSSatish Balay   PetscFunctionBegin;
859a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
860a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
861a7e14dcfSSatish Balay   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
862a7e14dcfSSatish Balay   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
863a7e14dcfSSatish Balay   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
864a7e14dcfSSatish Balay   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
865a7e14dcfSSatish Balay   nlsP->Diag = 0;
866a7e14dcfSSatish Balay   nlsP->M = 0;
867a7e14dcfSSatish Balay   PetscFunctionReturn(0);
868a7e14dcfSSatish Balay }
869a7e14dcfSSatish Balay 
870a7e14dcfSSatish Balay /*------------------------------------------------------------*/
871a7e14dcfSSatish Balay #undef __FUNCT__
872a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NLS"
873441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao)
874a7e14dcfSSatish Balay {
875a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
876a7e14dcfSSatish Balay   PetscErrorCode ierr;
877a7e14dcfSSatish Balay 
878a7e14dcfSSatish Balay   PetscFunctionBegin;
879a7e14dcfSSatish Balay   if (tao->setupcalled) {
880a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
881a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
882a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
883a7e14dcfSSatish Balay     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
884a7e14dcfSSatish Balay   }
885a7e14dcfSSatish Balay   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
886a7e14dcfSSatish Balay   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
887a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
888a7e14dcfSSatish Balay   PetscFunctionReturn(0);
889a7e14dcfSSatish Balay }
890a7e14dcfSSatish Balay 
891a7e14dcfSSatish Balay /*------------------------------------------------------------*/
892a7e14dcfSSatish Balay #undef __FUNCT__
893a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NLS"
8944416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptionItems *PetscOptionsObject,Tao tao)
895a7e14dcfSSatish Balay {
896a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
897a7e14dcfSSatish Balay   PetscErrorCode ierr;
898a7e14dcfSSatish Balay 
899a7e14dcfSSatish Balay   PetscFunctionBegin;
9001a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
901a7e14dcfSSatish 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);
902a7e14dcfSSatish 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);
903a7e14dcfSSatish 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);
904a7e14dcfSSatish 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);
905a7e14dcfSSatish 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);
90694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr);
90794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr);
90894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr);
90994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr);
91094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr);
91194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr);
91294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr);
91394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr);
91494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr);
91594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr);
91694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr);
91794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr);
91894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr);
91994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr);
92094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr);
92194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr);
92294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr);
92394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr);
92494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr);
92594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr);
92694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr);
92794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr);
92894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr);
92994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr);
93094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr);
93194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr);
93294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr);
93394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr);
93494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr);
93594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr);
93694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr);
93794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr);
93894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr);
93994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr);
94094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr);
94194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr);
94294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr);
94394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr);
94494ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr);
94594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr);
94694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr);
94794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr);
94894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr);
94994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr);
95094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr);
951a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
952a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
953a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
954a7e14dcfSSatish Balay   PetscFunctionReturn(0);
955a7e14dcfSSatish Balay }
956a7e14dcfSSatish Balay 
957a7e14dcfSSatish Balay 
958a7e14dcfSSatish Balay /*------------------------------------------------------------*/
959a7e14dcfSSatish Balay #undef __FUNCT__
960a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NLS"
961441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
962a7e14dcfSSatish Balay {
963a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
964a7e14dcfSSatish Balay   PetscInt       nrejects;
965a7e14dcfSSatish Balay   PetscBool      isascii;
966a7e14dcfSSatish Balay   PetscErrorCode ierr;
967a7e14dcfSSatish Balay 
968a7e14dcfSSatish Balay   PetscFunctionBegin;
969a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
970a7e14dcfSSatish Balay   if (isascii) {
971a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
972a7e14dcfSSatish Balay     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
973a7e14dcfSSatish Balay       ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr);
974a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
975a7e14dcfSSatish Balay     }
976a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
977a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
978a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
979a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
980a7e14dcfSSatish Balay 
981a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
982a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
983a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
984a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
985a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
986a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
987a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
988a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
989a7e14dcfSSatish Balay   }
990a7e14dcfSSatish Balay   PetscFunctionReturn(0);
991a7e14dcfSSatish Balay }
992a7e14dcfSSatish Balay 
993a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
9944aa34175SJason Sarich /*MC
9954aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
9964aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
9974aa34175SJason Sarich   system of equations to obtain the step diretion dk:
9984aa34175SJason Sarich               Hk dk = -gk
9994aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
10004aa34175SJason Sarich   solve
10014aa34175SJason Sarich               min_t f(xk + t d_k)
10024aa34175SJason Sarich 
10034aa34175SJason Sarich     Options Database Keys:
10044aa34175SJason Sarich + -tao_nls_ksp_type - "cg","nash","stcg","gltr","petsc"
10054aa34175SJason Sarich . -tao_nls_pc_type - "none","ahess","bfgs","petsc"
10064aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs"
10074aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation"
10084aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
10094aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
10104aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
10114aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
10124aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
10134aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
10144aa34175SJason Sarich . -tao_nls_pgfac - growth factor
10154aa34175SJason Sarich . -tao_nls_psfac - shrink factor
10164aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
10174aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
10184aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
10194aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
10204aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
10214aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
10224aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
10234aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
10244aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
10254aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
10264aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
10274aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
10284aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
10294aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
10304aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
10314aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
10324aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
10334aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
10344aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
10354aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
10364aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
10374aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
10384aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
10394aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
10401522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
10411522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
10421522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
10431522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
10441522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
10451522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
10461522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
10471eb8069cSJason Sarich 
10481eb8069cSJason Sarich   Level: beginner
10491522df2eSJason Sarich M*/
10504aa34175SJason Sarich 
1051a7e14dcfSSatish Balay #undef __FUNCT__
1052a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NLS"
1053728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
1054a7e14dcfSSatish Balay {
1055a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
10568caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
1057a7e14dcfSSatish Balay   PetscErrorCode ierr;
1058a7e14dcfSSatish Balay 
1059a7e14dcfSSatish Balay   PetscFunctionBegin;
10603c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1061a7e14dcfSSatish Balay 
1062a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
1063a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
1064a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
1065a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1066a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
1067a7e14dcfSSatish Balay 
10686552cf8aSJason Sarich   /* Override default settings (unless already changed) */
10696552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
10706552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
10716552cf8aSJason Sarich 
1072a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
1073a7e14dcfSSatish Balay 
1074a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
1075a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
1076a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
1077a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
1078a7e14dcfSSatish Balay 
1079a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
1080a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
1081a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
1082a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
1083a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
1084a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
1085a7e14dcfSSatish Balay 
1086a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
1087a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
1088a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
1089a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
1090a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
1091a7e14dcfSSatish Balay 
1092a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
1093a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
1094a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
1095a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
1096a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
1097a7e14dcfSSatish Balay 
1098a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
1099a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
1100a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
1101a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
1102a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
1103a7e14dcfSSatish Balay 
1104a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
1105a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
1106a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
1107a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
1108a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
1109a7e14dcfSSatish Balay 
1110a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
1111a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
1112a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
1113a7e14dcfSSatish Balay 
1114a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
1115a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
1116a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
1117a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
1118a7e14dcfSSatish Balay 
1119a7e14dcfSSatish Balay   nlsP->theta = 0.05;
1120a7e14dcfSSatish Balay 
1121a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
1122a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
1123a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
1124a7e14dcfSSatish Balay 
1125a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
1126a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
1127a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
1128a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
1129a7e14dcfSSatish Balay 
1130a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
1131a7e14dcfSSatish Balay 
1132a7e14dcfSSatish Balay   /*  Remaining parameters */
1133a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
1134a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
1135a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
1136a7e14dcfSSatish Balay 
1137a7e14dcfSSatish Balay   nlsP->ksp_type        = NLS_KSP_STCG;
1138a7e14dcfSSatish Balay   nlsP->pc_type         = NLS_PC_BFGS;
1139a7e14dcfSSatish Balay   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1140a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1141a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
1142a7e14dcfSSatish Balay 
1143a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1144a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1145441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
11465d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1147a7e14dcfSSatish Balay 
1148a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
1149a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
11505d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
1151a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1152a7e14dcfSSatish Balay }
1153a7e14dcfSSatish Balay 
1154a7e14dcfSSatish Balay #undef __FUNCT__
1155a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell"
1156a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
1157a7e14dcfSSatish Balay {
1158a7e14dcfSSatish Balay   PetscErrorCode ierr;
1159a7e14dcfSSatish Balay   Mat            M;
116087f595a5SBarry Smith 
1161a7e14dcfSSatish Balay   PetscFunctionBegin;
1162a7e14dcfSSatish Balay   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1163a7e14dcfSSatish Balay   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
1164a7e14dcfSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
1165a7e14dcfSSatish Balay   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
1166a7e14dcfSSatish Balay   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
1167a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1168a7e14dcfSSatish Balay }
1169