xref: /petsc/src/snes/impls/qn/qn.c (revision 88f769c5898b98f601d4274d5ebe4e042065520a)
10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
24b11644fSPeter Brune 
38e231d97SPeter Brune #define H(i,j)  qn->dXdFmat[i*qn->m + j]
48e231d97SPeter Brune 
56a6fc655SJed Brown const char *const SNESQNScaleTypes[] =        {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
66a6fc655SJed Brown const char *const SNESQNRestartTypes[] =      {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
7b8840d6bSPeter Brune const char *const SNESQNTypes[] =             {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0};
8b8840d6bSPeter Brune 
9b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS      = 0,
10b8840d6bSPeter Brune               SNES_QN_BROYDEN    = 1,
11b8840d6bSPeter Brune               SNES_QN_BADBROYDEN = 2} SNESQNType;
126bf1b2e5SPeter Brune 
134b11644fSPeter Brune typedef struct {
14b8840d6bSPeter Brune   Vec          *U;               /* Stored past states (vary from method to method) */
15b8840d6bSPeter Brune   Vec          *V;               /* Stored past states (vary from method to method) */
168e231d97SPeter Brune   PetscInt     m;                /* The number of kept previous steps */
178e231d97SPeter Brune   PetscScalar  *alpha, *beta;
188e231d97SPeter Brune   PetscScalar  *dXtdF, *dFtdX, *YtdX;
1918aa0c0cSPeter Brune   PetscBool    singlereduction;  /* Aggregated reduction implementation */
208e231d97SPeter Brune   PetscScalar  *dXdFmat;         /* A matrix of values for dX_i dot dF_j */
2144f7e39eSPeter Brune   PetscViewer  monitor;
226bf1b2e5SPeter Brune   PetscReal    powell_gamma;     /* Powell angle restart condition */
236bf1b2e5SPeter Brune   PetscReal    powell_downhill;  /* Powell descent restart condition */
24b21d5a53SPeter Brune   PetscReal    scaling;          /* scaling of H0 */
2588d374b2SPeter Brune 
26b8840d6bSPeter Brune   SNESQNType            type;             /* the type of quasi-newton method used */
27*88f769c5SPeter Brune   SNESQNScaleType       scale_type;       /* the type of scaling used */
280c777b0cSPeter Brune   SNESQNRestartType     restart_type;     /* determine the frequency and type of restart conditions */
299f83bee8SJed Brown } SNES_QN;
304b11644fSPeter Brune 
314b11644fSPeter Brune #undef __FUNCT__
32b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden"
33b8840d6bSPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold) {
344b11644fSPeter Brune 
354b11644fSPeter Brune   PetscErrorCode ierr;
364b11644fSPeter Brune 
379f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
384b11644fSPeter Brune 
39b8840d6bSPeter Brune   Vec W = snes->work[3];
400ecc9e1dSPeter Brune 
41b8840d6bSPeter Brune   Vec *U = qn->U;
42b8840d6bSPeter Brune   Vec *V = qn->V;
43b8840d6bSPeter Brune 
44b8840d6bSPeter Brune   /* ksp thing for jacobian scaling */
45b8840d6bSPeter Brune   KSPConvergedReason kspreason;
46b8840d6bSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
47b8840d6bSPeter Brune 
48b8840d6bSPeter Brune   PetscInt k,i,lits;
49b8840d6bSPeter Brune   PetscInt m = qn->m;
50b8840d6bSPeter Brune   PetscScalar gdot;
51b8840d6bSPeter Brune   PetscInt l = m;
52b8840d6bSPeter Brune 
53b8840d6bSPeter Brune   Mat jac,jac_pre;
54b8840d6bSPeter Brune   PetscFunctionBegin;
55b8840d6bSPeter Brune 
56b8840d6bSPeter Brune   if (it < m) l = it;
57b8840d6bSPeter Brune 
58b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
59b8840d6bSPeter Brune     ierr = SNESGetJacobian(snes,&jac,&jac_pre,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
60b8840d6bSPeter Brune     ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr);
61b8840d6bSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,D,W);CHKERRQ(ierr);
62b8840d6bSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
63b8840d6bSPeter Brune     if (kspreason < 0) {
64b8840d6bSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
65b8840d6bSPeter Brune         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
66b8840d6bSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
67b8840d6bSPeter Brune         PetscFunctionReturn(0);
68b8840d6bSPeter Brune       }
69b8840d6bSPeter Brune     }
70b8840d6bSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
71b8840d6bSPeter Brune     snes->linear_its += lits;
72b8840d6bSPeter Brune     ierr = VecCopy(W,Y);CHKERRQ(ierr);
73b8840d6bSPeter Brune   } else {
74b8840d6bSPeter Brune     ierr = VecCopy(D,Y);CHKERRQ(ierr);
75b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
76b8840d6bSPeter Brune   }
77b8840d6bSPeter Brune 
78b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
79b8840d6bSPeter Brune   if (it > 1) {
80b8840d6bSPeter Brune     for (i = 0; i < l-1; i++) {
81b8840d6bSPeter Brune       k = (it+i-l)%l;
82b8840d6bSPeter Brune       ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
83b8840d6bSPeter Brune       ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr);
84b8840d6bSPeter Brune       if (qn->monitor) {
85b8840d6bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
86b8840d6bSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
87b8840d6bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
88b8840d6bSPeter Brune       }
89b8840d6bSPeter Brune     }
90b8840d6bSPeter Brune   }
91b8840d6bSPeter Brune   if (it < m) l = it;
92b8840d6bSPeter Brune 
93b8840d6bSPeter Brune   /* set W to be the last step, post-linesearch */
94b8840d6bSPeter Brune   ierr = VecCopy(Xold,W);CHKERRQ(ierr);
95b8840d6bSPeter Brune   ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
96b8840d6bSPeter Brune 
97b8840d6bSPeter Brune   if (l > 0) {
98b8840d6bSPeter Brune     k = (it-1)%l;
99b8840d6bSPeter Brune     ierr = VecCopy(W,U[k]);CHKERRQ(ierr);
100b8840d6bSPeter Brune     ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr);
101b8840d6bSPeter Brune     ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr);
102b8840d6bSPeter Brune     ierr = VecCopy(Y,V[k]);CHKERRQ(ierr);
103b8840d6bSPeter Brune     ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr);
104b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
105b8840d6bSPeter Brune     ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr);
106b8840d6bSPeter Brune     if (qn->monitor) {
107b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
108b8840d6bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
109b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
110b8840d6bSPeter Brune     }
111b8840d6bSPeter Brune   }
112b8840d6bSPeter Brune   PetscFunctionReturn(0);
113b8840d6bSPeter Brune }
114b8840d6bSPeter Brune 
115b8840d6bSPeter Brune #undef __FUNCT__
116b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden"
117b8840d6bSPeter Brune PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) {
118b8840d6bSPeter Brune 
119b8840d6bSPeter Brune   PetscErrorCode ierr;
120b8840d6bSPeter Brune 
121b8840d6bSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
122b8840d6bSPeter Brune 
123b8840d6bSPeter Brune   Vec W = snes->work[3];
124b8840d6bSPeter Brune 
125b8840d6bSPeter Brune   Vec *U = qn->U;
126b8840d6bSPeter Brune   Vec *T = qn->V;
127b8840d6bSPeter Brune 
128b8840d6bSPeter Brune   /* ksp thing for jacobian scaling */
129b8840d6bSPeter Brune   KSPConvergedReason kspreason;
130b8840d6bSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
131b8840d6bSPeter Brune 
132b8840d6bSPeter Brune   PetscInt k, i, lits;
133b8840d6bSPeter Brune   PetscInt m = qn->m;
134b8840d6bSPeter Brune   PetscScalar gdot;
135b8840d6bSPeter Brune   PetscInt l = m;
136b8840d6bSPeter Brune 
137b8840d6bSPeter Brune   Mat jac, jac_pre;
138b8840d6bSPeter Brune   PetscFunctionBegin;
139b8840d6bSPeter Brune 
140b8840d6bSPeter Brune   if (it < m) l = it;
141b8840d6bSPeter Brune 
142b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
143b8840d6bSPeter Brune 
144b8840d6bSPeter Brune   if (l > 0) {
145b8840d6bSPeter Brune     k = (it-1)%l;
146b8840d6bSPeter Brune     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
147b8840d6bSPeter Brune     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
148b8840d6bSPeter Brune     ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr);
149b8840d6bSPeter Brune     ierr = VecCopy(D,T[k]);CHKERRQ(ierr);
150b8840d6bSPeter Brune     ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr);
151b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
152b8840d6bSPeter Brune     ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr);
153b8840d6bSPeter Brune     if (qn->monitor) {
154b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
155b8840d6bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
156b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
157b8840d6bSPeter Brune     }
158b8840d6bSPeter Brune   }
159b8840d6bSPeter Brune 
160b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
161b8840d6bSPeter Brune   if (it > 1) {
162b8840d6bSPeter Brune     for (i = 0; i < l-1; i++) {
163b8840d6bSPeter Brune       k = (it+i-l)%l;
164b8840d6bSPeter Brune       ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
165b8840d6bSPeter Brune       ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr);
166b8840d6bSPeter Brune       if (qn->monitor) {
167b8840d6bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
168b8840d6bSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
169b8840d6bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
170b8840d6bSPeter Brune       }
171b8840d6bSPeter Brune     }
172b8840d6bSPeter Brune   }
173b8840d6bSPeter Brune 
174b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
175b8840d6bSPeter Brune     ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
176b8840d6bSPeter Brune     ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr);
177b8840d6bSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr);
178b8840d6bSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
179b8840d6bSPeter Brune     if (kspreason < 0) {
180b8840d6bSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
181b8840d6bSPeter Brune         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
182b8840d6bSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
183b8840d6bSPeter Brune         PetscFunctionReturn(0);
184b8840d6bSPeter Brune       }
185b8840d6bSPeter Brune     }
186b8840d6bSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
187b8840d6bSPeter Brune     snes->linear_its += lits;
188b8840d6bSPeter Brune     ierr = VecCopy(W,Y);CHKERRQ(ierr);
189b8840d6bSPeter Brune   }  else {
190b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
191b8840d6bSPeter Brune   }
192b8840d6bSPeter Brune   PetscFunctionReturn(0);
193b8840d6bSPeter Brune }
194b8840d6bSPeter Brune 
195b8840d6bSPeter Brune #undef __FUNCT__
196b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS"
197b8840d6bSPeter Brune PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) {
198b8840d6bSPeter Brune 
199b8840d6bSPeter Brune   PetscErrorCode ierr;
200b8840d6bSPeter Brune 
201b8840d6bSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
202b8840d6bSPeter Brune 
203b8840d6bSPeter Brune   Vec W = snes->work[3];
204b8840d6bSPeter Brune 
205b8840d6bSPeter Brune   Vec *dX = qn->U;
206b8840d6bSPeter Brune   Vec *dF = qn->V;
2074b11644fSPeter Brune 
208bd052dfeSPeter Brune   PetscScalar *alpha    = qn->alpha;
209bd052dfeSPeter Brune   PetscScalar *beta     = qn->beta;
2108e231d97SPeter Brune   PetscScalar *dXtdF    = qn->dXtdF;
211b8840d6bSPeter Brune   PetscScalar *dFtdX    = qn->dFtdX;
2128e231d97SPeter Brune   PetscScalar *YtdX     = qn->YtdX;
213bd052dfeSPeter Brune 
2140ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
2150ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
2160ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
2170ecc9e1dSPeter Brune 
2188e231d97SPeter Brune   PetscInt k,i,j,g,lits;
2194b11644fSPeter Brune   PetscInt m = qn->m;
2204b11644fSPeter Brune   PetscScalar t;
2214b11644fSPeter Brune   PetscInt l = m;
2224b11644fSPeter Brune 
2238e231d97SPeter Brune   Mat jac,jac_pre;
2248e231d97SPeter Brune 
2254b11644fSPeter Brune   PetscFunctionBegin;
2264b11644fSPeter Brune 
2275ba6227bSPeter Brune   if (it < m) l = it;
2284b11644fSPeter Brune 
229b8840d6bSPeter Brune   if (it > 0) {
230b8840d6bSPeter Brune     k = (it - 1) % l;
231b8840d6bSPeter Brune     ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
232b8840d6bSPeter Brune     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
233b8840d6bSPeter Brune     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
234b8840d6bSPeter Brune     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
23518aa0c0cSPeter Brune     if (qn->singlereduction) {
236b8840d6bSPeter Brune       ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr);
237b8840d6bSPeter Brune       ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr);
238b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
239b8840d6bSPeter Brune         H(k, j) = dFtdX[j];
240b8840d6bSPeter Brune         H(j, k) = dXtdF[j];
241b8840d6bSPeter Brune       }
242b8840d6bSPeter Brune       /* copy back over to make the computation of alpha and beta easier */
243b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
244b8840d6bSPeter Brune         dXtdF[j] = H(j, j);
245b8840d6bSPeter Brune       }
246b8840d6bSPeter Brune     } else {
247b8840d6bSPeter Brune       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
248b8840d6bSPeter Brune     }
249b8840d6bSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_SHANNO) {
250b8840d6bSPeter Brune       PetscScalar dFtdF;
251b8840d6bSPeter Brune       ierr = VecDot(dF[k],dF[k],&dFtdF);CHKERRQ(ierr);
252b8840d6bSPeter Brune       qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(dFtdF);
253b8840d6bSPeter Brune     } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
254b8840d6bSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
255b8840d6bSPeter Brune     }
256b8840d6bSPeter Brune   }
257b8840d6bSPeter Brune 
258b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
259b8840d6bSPeter Brune 
260b8840d6bSPeter Brune   if (qn->singlereduction) {
261b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr);
2628e231d97SPeter Brune   }
2634b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
2644b11644fSPeter Brune   for (i=0;i<l;i++) {
265b21d5a53SPeter Brune     k = (it-i-1)%l;
26618aa0c0cSPeter Brune     if (qn->singlereduction) {
2678e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
2688e231d97SPeter Brune       t = YtdX[k];
2698e231d97SPeter Brune       for (j=0;j<i;j++) {
2708e231d97SPeter Brune         g = (it-j-1)%l;
2718e231d97SPeter Brune         t += -alpha[g]*H(g, k);
2728e231d97SPeter Brune       }
2738e231d97SPeter Brune       alpha[k] = t/H(k,k);
2748e231d97SPeter Brune     } else {
2753af51624SPeter Brune       ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
2768e231d97SPeter Brune       alpha[k] = t/dXtdF[k];
2778e231d97SPeter Brune     }
27844f7e39eSPeter Brune     if (qn->monitor) {
2793af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2805ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
2813af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
28244f7e39eSPeter Brune     }
2836bf1b2e5SPeter Brune     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
2844b11644fSPeter Brune   }
2854b11644fSPeter Brune 
2860c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
2878e231d97SPeter Brune     ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
2888e231d97SPeter Brune     ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr);
289b8840d6bSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr);
2900ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
2910ecc9e1dSPeter Brune     if (kspreason < 0) {
2920ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
2930ecc9e1dSPeter Brune         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
2940ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
2950ecc9e1dSPeter Brune         PetscFunctionReturn(0);
2960ecc9e1dSPeter Brune       }
2970ecc9e1dSPeter Brune     }
2980ecc9e1dSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2990ecc9e1dSPeter Brune     snes->linear_its += lits;
300b8840d6bSPeter Brune     ierr = VecCopy(W, Y);CHKERRQ(ierr);
3010ecc9e1dSPeter Brune   } else {
302b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
3030ecc9e1dSPeter Brune   }
30418aa0c0cSPeter Brune   if (qn->singlereduction) {
305b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
3068e231d97SPeter Brune   }
3074b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
308bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
309b21d5a53SPeter Brune     k = (it + i - l) % l;
31018aa0c0cSPeter Brune     if (qn->singlereduction) {
3118e231d97SPeter Brune       t = YtdX[k];
3128e231d97SPeter Brune       for (j = 0; j < i; j++) {
3138e231d97SPeter Brune         g = (it + j - l) % l;
3148e231d97SPeter Brune         t += (alpha[g] - beta[g])*H(k, g);
3158e231d97SPeter Brune       }
3168e231d97SPeter Brune       beta[k] = t / H(k, k);
3178e231d97SPeter Brune     } else {
3186bf1b2e5SPeter Brune       ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
3198e231d97SPeter Brune       beta[k] = t / dXtdF[k];
3208e231d97SPeter Brune     }
3213af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
32244f7e39eSPeter Brune     if (qn->monitor) {
3233af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3245ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
3253af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
32644f7e39eSPeter Brune     }
3274b11644fSPeter Brune   }
3284b11644fSPeter Brune   PetscFunctionReturn(0);
3294b11644fSPeter Brune }
3304b11644fSPeter Brune 
3314b11644fSPeter Brune #undef __FUNCT__
3324b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
3334b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
3344b11644fSPeter Brune {
3354b11644fSPeter Brune 
3364b11644fSPeter Brune   PetscErrorCode ierr;
3379f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
3384b11644fSPeter Brune 
33915f5eeeaSPeter Brune   Vec X,Xold;
340335fb975SPeter Brune   Vec F,B;
341335fb975SPeter Brune   Vec Y,FPC,D,Dold;
3423af51624SPeter Brune   SNESConvergedReason reason;
343b8840d6bSPeter Brune   PetscInt i, i_r;
3444b11644fSPeter Brune 
345335fb975SPeter Brune   PetscReal fnorm,xnorm,ynorm,gnorm;
3460c777b0cSPeter Brune   PetscBool lssucceed,powell,periodic;
3474b11644fSPeter Brune 
348b8840d6bSPeter Brune   PetscScalar DolddotD,DolddotDold,DdotD,YdotD;
3494b11644fSPeter Brune 
3500ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
3510ecc9e1dSPeter Brune 
3524b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
3534b11644fSPeter Brune   PetscFunctionBegin;
3544b11644fSPeter Brune 
3559f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
3563af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
3573af51624SPeter Brune   B             = snes->vec_rhs;
358b8840d6bSPeter Brune 
359b8840d6bSPeter Brune   X             = snes->vec_sol;        /* solution vector */
360335fb975SPeter Brune   Xold          = snes->work[0];
3613af51624SPeter Brune 
3623af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
363335fb975SPeter Brune   D             = snes->work[1];
364335fb975SPeter Brune   Dold          = snes->work[2];
3654b11644fSPeter Brune 
3664b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
3674b11644fSPeter Brune 
3684b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
3694b11644fSPeter Brune   snes->iter = 0;
3704b11644fSPeter Brune   snes->norm = 0.;
3714b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
372e4ed7901SPeter Brune   if (!snes->vec_func_init_set){
37315f5eeeaSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3744b11644fSPeter Brune     if (snes->domainerror) {
3754b11644fSPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
3764b11644fSPeter Brune       PetscFunctionReturn(0);
3774b11644fSPeter Brune     }
378e4ed7901SPeter Brune   } else {
379e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
380e4ed7901SPeter Brune   }
381e4ed7901SPeter Brune 
382e4ed7901SPeter Brune   if (!snes->norm_init_set) {
38315f5eeeaSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
3844b11644fSPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
385e4ed7901SPeter Brune   } else {
386e4ed7901SPeter Brune     fnorm = snes->norm_init;
387e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
388e4ed7901SPeter Brune   }
389e4ed7901SPeter Brune 
3904b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
3914b11644fSPeter Brune   snes->norm = fnorm;
3924b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
3934b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
3944b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
3954b11644fSPeter Brune 
3964b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
3974b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
3984b11644fSPeter Brune   /* test convergence */
3994b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4004b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
40170d3b23bSPeter Brune 
402*88f769c5SPeter Brune   /* composed solve */
403c40d0f55SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
404217b9c2eSPeter Brune     ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
405217b9c2eSPeter Brune     ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
40688d374b2SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
40788d374b2SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
40888d374b2SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
40988d374b2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
41088d374b2SPeter Brune       PetscFunctionReturn(0);
41188d374b2SPeter Brune     }
41288d374b2SPeter Brune     ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
41388d374b2SPeter Brune     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
41488d374b2SPeter Brune     ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
4156bf1b2e5SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
41688d374b2SPeter Brune   } else {
41788d374b2SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
41888d374b2SPeter Brune   }
41988d374b2SPeter Brune   ierr = VecCopy(Y, D);CHKERRQ(ierr);
4203af51624SPeter Brune 
421f8e15203SPeter Brune   /* scale the initial update */
4220c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
4230ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
4240ecc9e1dSPeter Brune   }
4250ecc9e1dSPeter Brune 
4265ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
427b8840d6bSPeter Brune     switch(qn->type) {
428b8840d6bSPeter Brune     case SNES_QN_BADBROYDEN:
429b8840d6bSPeter Brune       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
430b8840d6bSPeter Brune       break;
431b8840d6bSPeter Brune     case SNES_QN_BROYDEN:
432b8840d6bSPeter Brune       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
433b8840d6bSPeter Brune       break;
434b8840d6bSPeter Brune     case SNES_QN_LBFGS:
435b8840d6bSPeter Brune       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
436b8840d6bSPeter Brune       break;
437b8840d6bSPeter Brune     }
43870d3b23bSPeter Brune     /* line search for lambda */
43970d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
44088d374b2SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
44115f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
442f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
4439f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
4449f3a0142SPeter Brune     if (snes->domainerror) {
4459f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
4469f3a0142SPeter Brune       PetscFunctionReturn(0);
4479f3a0142SPeter Brune       }
448f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
4499f3a0142SPeter Brune     if (!lssucceed) {
4509f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
4519f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
4529f3a0142SPeter Brune         break;
4539f3a0142SPeter Brune       }
4549f3a0142SPeter Brune     }
455f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
4560c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
45705ee83feSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
4580ecc9e1dSPeter Brune     }
4593af51624SPeter Brune 
46088d374b2SPeter Brune     /* convergence monitoring */
461b21d5a53SPeter Brune     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
462b21d5a53SPeter Brune 
463360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
464360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
465360c497dSPeter Brune 
4668409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
4678409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
4684b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
469d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4704b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
4714b11644fSPeter Brune 
472c40d0f55SPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
473217b9c2eSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
474217b9c2eSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
47588d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
47688d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
47788d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
47888d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
47988d374b2SPeter Brune         PetscFunctionReturn(0);
48088d374b2SPeter Brune       }
48188d374b2SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
48288d374b2SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
48388d374b2SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
48488d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
48588d374b2SPeter Brune     } else {
48688d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
48788d374b2SPeter Brune     }
48888d374b2SPeter Brune 
4890c777b0cSPeter Brune     powell = PETSC_FALSE;
4900c777b0cSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
4916bf1b2e5SPeter Brune       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
4928e231d97SPeter Brune       ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
4938e231d97SPeter Brune       ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr);
4948e231d97SPeter Brune       ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr);
4958e231d97SPeter Brune       ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr);
4968e231d97SPeter Brune       ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
4978e231d97SPeter Brune       ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr);
4988e231d97SPeter Brune       ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr);
4998e231d97SPeter Brune       ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr);
5000c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
5010c777b0cSPeter Brune     }
5020c777b0cSPeter Brune     periodic = PETSC_FALSE;
503b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
504b8840d6bSPeter Brune       if (i_r>qn->m-1) periodic = PETSC_TRUE;
5050c777b0cSPeter Brune     }
5060c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
5070c777b0cSPeter Brune     if (powell || periodic) {
5085ba6227bSPeter Brune       if (qn->monitor) {
5095ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
510516fe3c3SPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr);
5115ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
5125ba6227bSPeter Brune       }
5135ba6227bSPeter Brune       i_r = -1;
5145ba6227bSPeter Brune       /* general purpose update */
5155ba6227bSPeter Brune       if (snes->ops->update) {
5165ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
5175ba6227bSPeter Brune       }
5180c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5190ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
5200ecc9e1dSPeter Brune       }
5210ecc9e1dSPeter Brune     }
52270d3b23bSPeter Brune     /* general purpose update */
52370d3b23bSPeter Brune     if (snes->ops->update) {
52470d3b23bSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
52570d3b23bSPeter Brune     }
5265ba6227bSPeter Brune   }
5274b11644fSPeter Brune   if (i == snes->max_its) {
5284b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
5294b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
5304b11644fSPeter Brune   }
5314b11644fSPeter Brune   PetscFunctionReturn(0);
5324b11644fSPeter Brune }
5334b11644fSPeter Brune 
5344b11644fSPeter Brune 
5354b11644fSPeter Brune #undef __FUNCT__
5364b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
5374b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
5384b11644fSPeter Brune {
5399f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
5404b11644fSPeter Brune   PetscErrorCode ierr;
541335fb975SPeter Brune 
5424b11644fSPeter Brune   PetscFunctionBegin;
543b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
544b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
5458e231d97SPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha,
5468e231d97SPeter Brune                       qn->m, PetscScalar, &qn->beta,
5478e231d97SPeter Brune                       qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr);
5488e231d97SPeter Brune 
54918aa0c0cSPeter Brune   if (qn->singlereduction) {
5508e231d97SPeter Brune     ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat,
5518e231d97SPeter Brune                         qn->m, PetscScalar, &qn->dFtdX,
5528e231d97SPeter Brune                         qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr);
5538e231d97SPeter Brune   }
554335fb975SPeter Brune   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
555335fb975SPeter Brune 
556335fb975SPeter Brune   /* set up the line search */
5570c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5588e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
5598e231d97SPeter Brune   }
5604b11644fSPeter Brune   PetscFunctionReturn(0);
5614b11644fSPeter Brune }
5624b11644fSPeter Brune 
5634b11644fSPeter Brune #undef __FUNCT__
5644b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
5654b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
5664b11644fSPeter Brune {
5674b11644fSPeter Brune   PetscErrorCode ierr;
5689f83bee8SJed Brown   SNES_QN *qn;
5694b11644fSPeter Brune   PetscFunctionBegin;
5704b11644fSPeter Brune   if (snes->data) {
5719f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
572b8840d6bSPeter Brune     if (qn->U) {
573b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
5744b11644fSPeter Brune     }
575b8840d6bSPeter Brune     if (qn->V) {
576b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
5774b11644fSPeter Brune     }
57818aa0c0cSPeter Brune     if (qn->singlereduction) {
5798e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
5808e231d97SPeter Brune     }
5818e231d97SPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr);
5824b11644fSPeter Brune   }
5834b11644fSPeter Brune   PetscFunctionReturn(0);
5844b11644fSPeter Brune }
5854b11644fSPeter Brune 
5864b11644fSPeter Brune #undef __FUNCT__
5874b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
5884b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
5894b11644fSPeter Brune {
5904b11644fSPeter Brune   PetscErrorCode ierr;
5914b11644fSPeter Brune   PetscFunctionBegin;
5924b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
5934b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
5949c05af5eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
5954b11644fSPeter Brune   PetscFunctionReturn(0);
5964b11644fSPeter Brune }
5974b11644fSPeter Brune 
5984b11644fSPeter Brune #undef __FUNCT__
5994b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
6004b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
6014b11644fSPeter Brune {
6024b11644fSPeter Brune 
6034b11644fSPeter Brune   PetscErrorCode ierr;
6049f83bee8SJed Brown   SNES_QN    *qn;
605*88f769c5SPeter Brune   PetscBool  monflg = PETSC_FALSE,flg;
606f1c6b773SPeter Brune   SNESLineSearch linesearch;
607*88f769c5SPeter Brune   SNESQNRestartType rtype;
608*88f769c5SPeter Brune   SNESQNScaleType   stype;
6094b11644fSPeter Brune   PetscFunctionBegin;
6104b11644fSPeter Brune 
6119f83bee8SJed Brown   qn    = (SNES_QN*)snes->data;
612*88f769c5SPeter Brune   rtype = qn->restart_type;
613*88f769c5SPeter Brune   stype = qn->scale_type;
6144b11644fSPeter Brune 
6154b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
6165b4627e8SPeter Brune   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr);
6175b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
6185b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
6195b4627e8SPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
6204d116c7dSPeter Brune   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr);
621*88f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
622*88f769c5SPeter Brune   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
623*88f769c5SPeter Brune 
624*88f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
625*88f769c5SPeter Brune   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
626*88f769c5SPeter Brune 
627b8840d6bSPeter Brune   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,
628b8840d6bSPeter Brune                           (PetscEnum)qn->type,(PetscEnum*)&qn->type,PETSC_NULL);CHKERRQ(ierr);
6294b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
6309e764e56SPeter Brune   if (!snes->linesearch) {
631f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
6321a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
6339e764e56SPeter Brune   }
63444f7e39eSPeter Brune   if (monflg) {
63544f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
63644f7e39eSPeter Brune   }
6374b11644fSPeter Brune   PetscFunctionReturn(0);
6384b11644fSPeter Brune }
6394b11644fSPeter Brune 
64092c02d66SPeter Brune 
6410c777b0cSPeter Brune #undef __FUNCT__
6420c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType"
6430c777b0cSPeter Brune /*@
6440c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
6450c777b0cSPeter Brune 
6460c777b0cSPeter Brune     Logically Collective on SNES
6470c777b0cSPeter Brune 
6480c777b0cSPeter Brune     Input Parameters:
6490c777b0cSPeter Brune +   snes - the iterative context
6500c777b0cSPeter Brune -   rtype - restart type
6510c777b0cSPeter Brune 
6520c777b0cSPeter Brune     Options Database:
6530c777b0cSPeter Brune +   -snes_qn_restart_type<powell,periodic,none> - set the restart type
654b8840d6bSPeter Brune -   -snes_qn_restart[10] - sets the number of iterations before restart for periodic
6550c777b0cSPeter Brune 
6560c777b0cSPeter Brune     Level: intermediate
6570c777b0cSPeter Brune 
6580c777b0cSPeter Brune     SNESQNRestartTypes:
6590c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
6600c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
6610c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
6620c777b0cSPeter Brune 
6630c777b0cSPeter Brune     Notes:
6640c777b0cSPeter Brune     The default line search used is the L2 line search and it requires two additional function evaluations.
6650c777b0cSPeter Brune 
6660c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch
6670c777b0cSPeter Brune @*/
6680c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) {
6690c777b0cSPeter Brune   PetscErrorCode ierr;
6700c777b0cSPeter Brune   PetscFunctionBegin;
6710c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6720c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
6730c777b0cSPeter Brune   PetscFunctionReturn(0);
6740c777b0cSPeter Brune }
6750c777b0cSPeter Brune 
6760c777b0cSPeter Brune #undef __FUNCT__
6770c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType"
6780c777b0cSPeter Brune /*@
6790c777b0cSPeter Brune     SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN.
6800c777b0cSPeter Brune 
6810c777b0cSPeter Brune     Logically Collective on SNES
6820c777b0cSPeter Brune 
6830c777b0cSPeter Brune     Input Parameters:
6840c777b0cSPeter Brune +   snes - the iterative context
6850c777b0cSPeter Brune -   stype - scale type
6860c777b0cSPeter Brune 
6870c777b0cSPeter Brune     Options Database:
6880c777b0cSPeter Brune .   -snes_qn_scale_type<shanno,none,linesearch,jacobian>
6890c777b0cSPeter Brune 
6900c777b0cSPeter Brune     Level: intermediate
6910c777b0cSPeter Brune 
6920c777b0cSPeter Brune     SNESQNSelectTypes:
6930c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
6940c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
6950c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
6960c777b0cSPeter Brune -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
6970c777b0cSPeter Brune 
6980c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
6990c777b0cSPeter Brune @*/
7000c777b0cSPeter Brune 
7010c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) {
7020c777b0cSPeter Brune   PetscErrorCode ierr;
7030c777b0cSPeter Brune   PetscFunctionBegin;
7040c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7050c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
7060c777b0cSPeter Brune   PetscFunctionReturn(0);
7070c777b0cSPeter Brune }
7080c777b0cSPeter Brune 
7090c777b0cSPeter Brune EXTERN_C_BEGIN
7100c777b0cSPeter Brune #undef __FUNCT__
7110c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN"
7120c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) {
7130c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
7140c777b0cSPeter Brune   PetscFunctionBegin;
7150c777b0cSPeter Brune   qn->scale_type = stype;
7160c777b0cSPeter Brune   PetscFunctionReturn(0);
7170c777b0cSPeter Brune }
7180c777b0cSPeter Brune 
7190c777b0cSPeter Brune #undef __FUNCT__
7200c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN"
7210c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) {
7220c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
7230c777b0cSPeter Brune   PetscFunctionBegin;
7240c777b0cSPeter Brune   qn->restart_type = rtype;
7250c777b0cSPeter Brune   PetscFunctionReturn(0);
7260c777b0cSPeter Brune }
7270c777b0cSPeter Brune EXTERN_C_END
7280c777b0cSPeter Brune 
7294b11644fSPeter Brune /* -------------------------------------------------------------------------- */
7304b11644fSPeter Brune /*MC
7314b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
7324b11644fSPeter Brune 
7336cc8130cSPeter Brune       Options Database:
7346cc8130cSPeter Brune 
7356cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
7361867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
7371867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
73872835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
73944f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
7406cc8130cSPeter Brune 
741b8840d6bSPeter Brune       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
742b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
743b8840d6bSPeter Brune       updates.
7446cc8130cSPeter Brune 
7451867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
7461867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
7471867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
7481867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
7491867fe5bSPeter Brune 
7506cc8130cSPeter Brune       References:
7516cc8130cSPeter Brune 
7526cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
7536cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
7546cc8130cSPeter Brune 
755b8840d6bSPeter Brune       Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker
756b8840d6bSPeter Brune       SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
757b8840d6bSPeter Brune 
7584b11644fSPeter Brune 
7594b11644fSPeter Brune       Level: beginner
7604b11644fSPeter Brune 
7614b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
7626cc8130cSPeter Brune 
7634b11644fSPeter Brune M*/
7644b11644fSPeter Brune EXTERN_C_BEGIN
7654b11644fSPeter Brune #undef __FUNCT__
7664b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
7674b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
7684b11644fSPeter Brune {
7694b11644fSPeter Brune 
7704b11644fSPeter Brune   PetscErrorCode ierr;
7719f83bee8SJed Brown   SNES_QN *qn;
7724b11644fSPeter Brune 
7734b11644fSPeter Brune   PetscFunctionBegin;
7744b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
7754b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
7764b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
7774b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
7784b11644fSPeter Brune   snes->ops->view            = 0;
7794b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
7804b11644fSPeter Brune 
78142f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
78242f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
78342f4f86dSBarry Smith 
78488976e71SPeter Brune   if (!snes->tolerancesset) {
7850e444f03SPeter Brune     snes->max_funcs = 30000;
7860e444f03SPeter Brune     snes->max_its   = 10000;
78788976e71SPeter Brune   }
7880e444f03SPeter Brune 
7899f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
7904b11644fSPeter Brune   snes->data = (void *) qn;
7910ecc9e1dSPeter Brune   qn->m               = 10;
792b21d5a53SPeter Brune   qn->scaling         = 1.0;
793b8840d6bSPeter Brune   qn->U               = PETSC_NULL;
794b8840d6bSPeter Brune   qn->V               = PETSC_NULL;
7958e231d97SPeter Brune   qn->dXtdF           = PETSC_NULL;
7968e231d97SPeter Brune   qn->dFtdX           = PETSC_NULL;
7978e231d97SPeter Brune   qn->dXdFmat         = PETSC_NULL;
79844f7e39eSPeter Brune   qn->monitor         = PETSC_NULL;
79918aa0c0cSPeter Brune   qn->singlereduction = PETSC_FALSE;
800b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
8010c777b0cSPeter Brune   qn->scale_type      = SNES_QN_SCALE_SHANNO;
8020c777b0cSPeter Brune   qn->restart_type    = SNES_QN_RESTART_POWELL;
803b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
804ea630c6eSPeter Brune 
805*88f769c5SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetScaleType_C","SNESQNSetScaleType_QN",SNESQNSetScaleType_QN);CHKERRQ(ierr);
806*88f769c5SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetRestartType_C","SNESQNSetRestartType_QN",SNESQNSetRestartType_QN);CHKERRQ(ierr);
807*88f769c5SPeter Brune 
8084b11644fSPeter Brune   PetscFunctionReturn(0);
8094b11644fSPeter Brune }
8108e231d97SPeter Brune 
8114b11644fSPeter Brune EXTERN_C_END
812