xref: /petsc/src/snes/impls/qn/qn.c (revision ddd40ce523a7cea6c5531ba94c3449af29fb68c1)
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,
110ad7597dSKarl Rupp               SNES_QN_BADBROYDEN = 2
120ad7597dSKarl Rupp              } SNESQNType;
136bf1b2e5SPeter Brune 
144b11644fSPeter Brune typedef struct {
15b8840d6bSPeter Brune   Vec               *U;                   /* Stored past states (vary from method to method) */
16b8840d6bSPeter Brune   Vec               *V;                   /* Stored past states (vary from method to method) */
178e231d97SPeter Brune   PetscInt          m;                    /* The number of kept previous steps */
188e231d97SPeter Brune   PetscScalar       *alpha, *beta;
198e231d97SPeter Brune   PetscScalar       *dXtdF, *dFtdX, *YtdX;
2018aa0c0cSPeter Brune   PetscBool         singlereduction;      /* Aggregated reduction implementation */
218e231d97SPeter Brune   PetscScalar       *dXdFmat;             /* A matrix of values for dX_i dot dF_j */
2244f7e39eSPeter Brune   PetscViewer       monitor;
236bf1b2e5SPeter Brune   PetscReal         powell_gamma;         /* Powell angle restart condition */
246bf1b2e5SPeter Brune   PetscReal         powell_downhill;      /* Powell descent restart condition */
25b21d5a53SPeter Brune   PetscReal         scaling;              /* scaling of H0 */
26b8840d6bSPeter Brune   SNESQNType        type;                 /* the type of quasi-newton method used */
2788f769c5SPeter 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"
332150357eSBarry Smith PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold)
342150357eSBarry Smith {
354b11644fSPeter Brune   PetscErrorCode     ierr;
369f83bee8SJed Brown   SNES_QN            *qn = (SNES_QN*)snes->data;
37b8840d6bSPeter Brune   Vec                W   = snes->work[3];
38b8840d6bSPeter Brune   Vec                *U  = qn->U;
39b8840d6bSPeter Brune   Vec                *V  = qn->V;
40b8840d6bSPeter Brune   KSPConvergedReason kspreason;
41b8840d6bSPeter Brune   PetscInt           k,i,lits;
42b8840d6bSPeter Brune   PetscInt           m = qn->m;
43b8840d6bSPeter Brune   PetscScalar        gdot;
44b8840d6bSPeter Brune   PetscInt           l = m;
452150357eSBarry Smith 
46b8840d6bSPeter Brune   PetscFunctionBegin;
47b8840d6bSPeter Brune   if (it < m) l = it;
48b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
49d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr);
50b8840d6bSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
51b8840d6bSPeter Brune     if (kspreason < 0) {
52b8840d6bSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
53b8840d6bSPeter 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);
54b8840d6bSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
55b8840d6bSPeter Brune         PetscFunctionReturn(0);
56b8840d6bSPeter Brune       }
57b8840d6bSPeter Brune     }
58b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
59b8840d6bSPeter Brune     snes->linear_its += lits;
60b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
61b8840d6bSPeter Brune   } else {
62b8840d6bSPeter Brune     ierr = VecCopy(D,Y);CHKERRQ(ierr);
63b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
64b8840d6bSPeter Brune   }
65b8840d6bSPeter Brune 
66b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
67b8840d6bSPeter Brune   if (it > 1) {
68b8840d6bSPeter Brune     for (i = 0; i < l-1; i++) {
69b8840d6bSPeter Brune       k    = (it+i-l)%l;
70b8840d6bSPeter Brune       ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
71b8840d6bSPeter Brune       ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr);
72b8840d6bSPeter Brune       if (qn->monitor) {
73b8840d6bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
74b8840d6bSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
75b8840d6bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
76b8840d6bSPeter Brune       }
77b8840d6bSPeter Brune     }
78b8840d6bSPeter Brune   }
79b8840d6bSPeter Brune   if (it < m) l = it;
80b8840d6bSPeter Brune 
81b8840d6bSPeter Brune   /* set W to be the last step, post-linesearch */
82b8840d6bSPeter Brune   ierr = VecCopy(Xold,W);CHKERRQ(ierr);
83b8840d6bSPeter Brune   ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
84b8840d6bSPeter Brune 
85b8840d6bSPeter Brune   if (l > 0) {
86b8840d6bSPeter Brune     k    = (it-1)%l;
87b8840d6bSPeter Brune     ierr = VecCopy(W,U[k]);CHKERRQ(ierr);
88b8840d6bSPeter Brune     ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr);
89b8840d6bSPeter Brune     ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr);
90b8840d6bSPeter Brune     ierr = VecCopy(Y,V[k]);CHKERRQ(ierr);
91b8840d6bSPeter Brune     ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr);
92b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
93b8840d6bSPeter Brune     ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr);
94b8840d6bSPeter Brune     if (qn->monitor) {
95b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
96b8840d6bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
97b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
98b8840d6bSPeter Brune     }
99b8840d6bSPeter Brune   }
100b8840d6bSPeter Brune   PetscFunctionReturn(0);
101b8840d6bSPeter Brune }
102b8840d6bSPeter Brune 
103b8840d6bSPeter Brune #undef __FUNCT__
104b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden"
1050adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1060adebc6cSBarry Smith {
107b8840d6bSPeter Brune   PetscErrorCode ierr;
108b8840d6bSPeter Brune   SNES_QN        *qn = (SNES_QN*)snes->data;
109b8840d6bSPeter Brune   Vec            W   = snes->work[3];
110b8840d6bSPeter Brune   Vec            *U  = qn->U;
111b8840d6bSPeter Brune   Vec            *T  = qn->V;
112b8840d6bSPeter Brune 
113b8840d6bSPeter Brune   /* ksp thing for jacobian scaling */
114b8840d6bSPeter Brune   KSPConvergedReason kspreason;
115b8840d6bSPeter Brune   PetscInt           k, i, lits;
116b8840d6bSPeter Brune   PetscInt           m = qn->m;
117b8840d6bSPeter Brune   PetscScalar        gdot;
118b8840d6bSPeter Brune   PetscInt           l = m;
1190adebc6cSBarry Smith 
120b8840d6bSPeter Brune   PetscFunctionBegin;
121b8840d6bSPeter Brune   if (it < m) l = it;
122b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
123b8840d6bSPeter Brune   if (l > 0) {
124b8840d6bSPeter Brune     k    = (it-1)%l;
125b8840d6bSPeter Brune     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
126b8840d6bSPeter Brune     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
127b8840d6bSPeter Brune     ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr);
128b8840d6bSPeter Brune     ierr = VecCopy(D,T[k]);CHKERRQ(ierr);
129b8840d6bSPeter Brune     ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr);
130b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
131b8840d6bSPeter Brune     ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr);
132b8840d6bSPeter Brune     if (qn->monitor) {
133b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
134b8840d6bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
135b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
136b8840d6bSPeter Brune     }
137b8840d6bSPeter Brune   }
138b8840d6bSPeter Brune 
139b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
140b8840d6bSPeter Brune   if (it > 1) {
141b8840d6bSPeter Brune     for (i = 0; i < l-1; i++) {
142b8840d6bSPeter Brune       k    = (it+i-l)%l;
143b8840d6bSPeter Brune       ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
144b8840d6bSPeter Brune       ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr);
145b8840d6bSPeter Brune       if (qn->monitor) {
146b8840d6bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
147b8840d6bSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
148b8840d6bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
149b8840d6bSPeter Brune       }
150b8840d6bSPeter Brune     }
151b8840d6bSPeter Brune   }
152b8840d6bSPeter Brune 
153b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
154d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
155b8840d6bSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
156b8840d6bSPeter Brune     if (kspreason < 0) {
157b8840d6bSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
158b8840d6bSPeter 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);
159b8840d6bSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
160b8840d6bSPeter Brune         PetscFunctionReturn(0);
161b8840d6bSPeter Brune       }
162b8840d6bSPeter Brune     }
163b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
164b8840d6bSPeter Brune     snes->linear_its += lits;
165b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
166b8840d6bSPeter Brune   } else {
167b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
168b8840d6bSPeter Brune   }
169b8840d6bSPeter Brune   PetscFunctionReturn(0);
170b8840d6bSPeter Brune }
171b8840d6bSPeter Brune 
172b8840d6bSPeter Brune #undef __FUNCT__
173b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS"
1740adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1750adebc6cSBarry Smith {
176b8840d6bSPeter Brune   PetscErrorCode ierr;
177b8840d6bSPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
178b8840d6bSPeter Brune   Vec            W      = snes->work[3];
179b8840d6bSPeter Brune   Vec            *dX    = qn->U;
180b8840d6bSPeter Brune   Vec            *dF    = qn->V;
181bd052dfeSPeter Brune   PetscScalar    *alpha = qn->alpha;
182bd052dfeSPeter Brune   PetscScalar    *beta  = qn->beta;
1838e231d97SPeter Brune   PetscScalar    *dXtdF = qn->dXtdF;
184b8840d6bSPeter Brune   PetscScalar    *dFtdX = qn->dFtdX;
1858e231d97SPeter Brune   PetscScalar    *YtdX  = qn->YtdX;
186bd052dfeSPeter Brune 
1870ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
1880ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
1898e231d97SPeter Brune   PetscInt           k,i,j,g,lits;
1904b11644fSPeter Brune   PetscInt           m = qn->m;
1914b11644fSPeter Brune   PetscScalar        t;
1924b11644fSPeter Brune   PetscInt           l = m;
1938e231d97SPeter Brune 
1944b11644fSPeter Brune   PetscFunctionBegin;
1955ba6227bSPeter Brune   if (it < m) l = it;
196b8840d6bSPeter Brune   if (it > 0) {
197b8840d6bSPeter Brune     k    = (it - 1) % l;
198b8840d6bSPeter Brune     ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
199b8840d6bSPeter Brune     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
200b8840d6bSPeter Brune     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
201b8840d6bSPeter Brune     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
20218aa0c0cSPeter Brune     if (qn->singlereduction) {
203b8840d6bSPeter Brune       ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr);
204b8840d6bSPeter Brune       ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr);
205b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
206b8840d6bSPeter Brune         H(k, j) = dFtdX[j];
207b8840d6bSPeter Brune         H(j, k) = dXtdF[j];
208b8840d6bSPeter Brune       }
209b8840d6bSPeter Brune       /* copy back over to make the computation of alpha and beta easier */
2101aa26658SKarl Rupp       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
211b8840d6bSPeter Brune     } else {
212b8840d6bSPeter Brune       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
213b8840d6bSPeter Brune     }
214b8840d6bSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_SHANNO) {
21567392de3SBarry Smith       PetscReal dFtdF;
21667392de3SBarry Smith       ierr        = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr);
21767392de3SBarry Smith       qn->scaling = PetscRealPart(dXtdF[k])/dFtdF;
218b8840d6bSPeter Brune     } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
219b8840d6bSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
220b8840d6bSPeter Brune     }
221b8840d6bSPeter Brune   }
222b8840d6bSPeter Brune 
223b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
224b8840d6bSPeter Brune 
225b8840d6bSPeter Brune   if (qn->singlereduction) {
226b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr);
2278e231d97SPeter Brune   }
2284b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
2294b11644fSPeter Brune   for (i=0; i<l; i++) {
230b21d5a53SPeter Brune     k = (it-i-1)%l;
23118aa0c0cSPeter Brune     if (qn->singlereduction) {
2328e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
2338e231d97SPeter Brune       t = YtdX[k];
2348e231d97SPeter Brune       for (j=0; j<i; j++) {
2358e231d97SPeter Brune         g  = (it-j-1)%l;
2368e231d97SPeter Brune         t += -alpha[g]*H(g, k);
2378e231d97SPeter Brune       }
2388e231d97SPeter Brune       alpha[k] = t/H(k,k);
2398e231d97SPeter Brune     } else {
2403af51624SPeter Brune       ierr     = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
2418e231d97SPeter Brune       alpha[k] = t/dXtdF[k];
2428e231d97SPeter Brune     }
24344f7e39eSPeter Brune     if (qn->monitor) {
2443af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2455ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
2463af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
24744f7e39eSPeter Brune     }
2486bf1b2e5SPeter Brune     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
2494b11644fSPeter Brune   }
2504b11644fSPeter Brune 
2510c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
252d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
2530ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
2540ecc9e1dSPeter Brune     if (kspreason < 0) {
2550ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
2560ecc9e1dSPeter 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);
2570ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
2580ecc9e1dSPeter Brune         PetscFunctionReturn(0);
2590ecc9e1dSPeter Brune       }
2600ecc9e1dSPeter Brune     }
2610ecc9e1dSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2620ecc9e1dSPeter Brune     snes->linear_its += lits;
263b8840d6bSPeter Brune     ierr              = VecCopy(W, Y);CHKERRQ(ierr);
2640ecc9e1dSPeter Brune   } else {
265b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
2660ecc9e1dSPeter Brune   }
26718aa0c0cSPeter Brune   if (qn->singlereduction) {
268b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
2698e231d97SPeter Brune   }
2704b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
271bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
272b21d5a53SPeter Brune     k = (it + i - l) % l;
27318aa0c0cSPeter Brune     if (qn->singlereduction) {
2748e231d97SPeter Brune       t = YtdX[k];
2758e231d97SPeter Brune       for (j = 0; j < i; j++) {
2768e231d97SPeter Brune         g  = (it + j - l) % l;
2778e231d97SPeter Brune         t += (alpha[g] - beta[g])*H(k, g);
2788e231d97SPeter Brune       }
2798e231d97SPeter Brune       beta[k] = t / H(k, k);
2808e231d97SPeter Brune     } else {
2816bf1b2e5SPeter Brune       ierr    = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
2828e231d97SPeter Brune       beta[k] = t / dXtdF[k];
2838e231d97SPeter Brune     }
28422d28d08SBarry Smith     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr);
28544f7e39eSPeter Brune     if (qn->monitor) {
2863af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2875ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
2883af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
28944f7e39eSPeter Brune     }
2904b11644fSPeter Brune   }
2914b11644fSPeter Brune   PetscFunctionReturn(0);
2924b11644fSPeter Brune }
2934b11644fSPeter Brune 
2944b11644fSPeter Brune #undef __FUNCT__
2954b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
2964b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
2974b11644fSPeter Brune {
2984b11644fSPeter Brune   PetscErrorCode      ierr;
2999f83bee8SJed Brown   SNES_QN             *qn = (SNES_QN*) snes->data;
30015f5eeeaSPeter Brune   Vec                 X,Xold;
30147f26062SPeter Brune   Vec                 F;
30247f26062SPeter Brune   Vec                 Y,D,Dold;
303b8840d6bSPeter Brune   PetscInt            i, i_r;
304335fb975SPeter Brune   PetscReal           fnorm,xnorm,ynorm,gnorm;
3050c777b0cSPeter Brune   PetscBool           lssucceed,powell,periodic;
306b8840d6bSPeter Brune   PetscScalar         DolddotD,DolddotDold,DdotD,YdotD;
3070ecc9e1dSPeter Brune   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
308b7281c8aSPeter Brune   SNESConvergedReason reason;
3090ecc9e1dSPeter Brune 
3104b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
3114b11644fSPeter Brune 
3126e111a19SKarl Rupp   PetscFunctionBegin;
3139f3a0142SPeter Brune   F = snes->vec_func;                   /* residual vector */
3143af51624SPeter Brune   Y = snes->vec_sol_update;             /* search direction generated by J^-1D*/
315b8840d6bSPeter Brune 
316b8840d6bSPeter Brune   X    = snes->vec_sol;                 /* solution vector */
317335fb975SPeter Brune   Xold = snes->work[0];
3183af51624SPeter Brune 
3193af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
320335fb975SPeter Brune   D    = snes->work[1];
321335fb975SPeter Brune   Dold = snes->work[2];
3224b11644fSPeter Brune 
3234b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
3244b11644fSPeter Brune 
325ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3264b11644fSPeter Brune   snes->iter = 0;
3274b11644fSPeter Brune   snes->norm = 0.;
328ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
32947f26062SPeter Brune 
330b28a06ddSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
33147f26062SPeter Brune     ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr);
332b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
333b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
334b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
335b7281c8aSPeter Brune       PetscFunctionReturn(0);
336b7281c8aSPeter Brune     }
33747f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
33847f26062SPeter Brune   } else {
339e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
34015f5eeeaSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3414b11644fSPeter Brune       if (snes->domainerror) {
3424b11644fSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
3434b11644fSPeter Brune         PetscFunctionReturn(0);
3444b11644fSPeter Brune       }
3451aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
346e4ed7901SPeter Brune     if (!snes->norm_init_set) {
34747f26062SPeter Brune       /* convergence test */
34847f26062SPeter Brune       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
349189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
350189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
351189a9710SBarry Smith         PetscFunctionReturn(0);
352189a9710SBarry Smith       }
353e4ed7901SPeter Brune     } else {
354e4ed7901SPeter Brune       fnorm               = snes->norm_init;
355e4ed7901SPeter Brune       snes->norm_init_set = PETSC_FALSE;
356e4ed7901SPeter Brune     }
35747f26062SPeter Brune   }
358b28a06ddSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
35947f26062SPeter Brune       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
360b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
361b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
362b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
363b7281c8aSPeter Brune         PetscFunctionReturn(0);
364b7281c8aSPeter Brune       }
36547f26062SPeter Brune   } else {
36647f26062SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
36747f26062SPeter Brune   }
368b28a06ddSPeter Brune 
369ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3704b11644fSPeter Brune   snes->norm = fnorm;
371ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
372a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
3734b11644fSPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
3744b11644fSPeter Brune 
3754b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
3764b11644fSPeter Brune   snes->ttol = fnorm*snes->rtol;
3774b11644fSPeter Brune   /* test convergence */
3784b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3794b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
38070d3b23bSPeter Brune 
3813cf07b75SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
3823cf07b75SPeter Brune     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
3833cf07b75SPeter Brune     ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
3843cf07b75SPeter Brune     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
385*ddd40ce5SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
386*ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
387*ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
388*ddd40ce5SPeter Brune       PetscFunctionReturn(0);
389*ddd40ce5SPeter Brune     }
390*ddd40ce5SPeter Brune     ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
3913cf07b75SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
3923cf07b75SPeter Brune   }
3933cf07b75SPeter Brune 
394f8e15203SPeter Brune   /* scale the initial update */
3950c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
3960ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
3978cba9625SJed Brown     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
3980ecc9e1dSPeter Brune   }
3990ecc9e1dSPeter Brune 
4005ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
401b8840d6bSPeter Brune     switch (qn->type) {
402b8840d6bSPeter Brune     case SNES_QN_BADBROYDEN:
403b8840d6bSPeter Brune       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
404b8840d6bSPeter Brune       break;
405b8840d6bSPeter Brune     case SNES_QN_BROYDEN:
406b8840d6bSPeter Brune       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
407b8840d6bSPeter Brune       break;
408b8840d6bSPeter Brune     case SNES_QN_LBFGS:
409b8840d6bSPeter Brune       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
410b8840d6bSPeter Brune       break;
411b8840d6bSPeter Brune     }
41270d3b23bSPeter Brune     /* line search for lambda */
41370d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
41488d374b2SPeter Brune     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
41515f5eeeaSPeter Brune     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
416f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
4179f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
4189f3a0142SPeter Brune     if (snes->domainerror) {
4199f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
4209f3a0142SPeter Brune       PetscFunctionReturn(0);
4219f3a0142SPeter Brune     }
422f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
4239f3a0142SPeter Brune     if (!lssucceed) {
4249f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
4259f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
4269f3a0142SPeter Brune         break;
4279f3a0142SPeter Brune       }
4289f3a0142SPeter Brune     }
429f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
4300c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
43105ee83feSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
4320ecc9e1dSPeter Brune     }
4333af51624SPeter Brune 
43488d374b2SPeter Brune     /* convergence monitoring */
435b21d5a53SPeter 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);
436b21d5a53SPeter Brune 
437b28a06ddSPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
438b28a06ddSPeter Brune       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
439b28a06ddSPeter Brune       ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
440b28a06ddSPeter Brune       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
441*ddd40ce5SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
442*ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
443*ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
444*ddd40ce5SPeter Brune         PetscFunctionReturn(0);
445*ddd40ce5SPeter Brune       }
446*ddd40ce5SPeter Brune       ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
447b28a06ddSPeter Brune     }
448b28a06ddSPeter Brune 
449360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
450360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
451360c497dSPeter Brune 
452a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
4538409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
4544b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
455d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4564b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
457b28a06ddSPeter Brune     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
45847f26062SPeter Brune       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
459b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
460b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
461b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
462b7281c8aSPeter Brune         PetscFunctionReturn(0);
463b7281c8aSPeter Brune       }
46488d374b2SPeter Brune     } else {
46588d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
46688d374b2SPeter Brune     }
46788d374b2SPeter Brune 
4680c777b0cSPeter Brune     powell = PETSC_FALSE;
4690c777b0cSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
4706bf1b2e5SPeter Brune       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
4718e231d97SPeter Brune       ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
4728e231d97SPeter Brune       ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr);
4738e231d97SPeter Brune       ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr);
4748e231d97SPeter Brune       ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr);
4758e231d97SPeter Brune       ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
4768e231d97SPeter Brune       ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr);
4778e231d97SPeter Brune       ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr);
4788e231d97SPeter Brune       ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr);
4790c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
4800c777b0cSPeter Brune     }
4810c777b0cSPeter Brune     periodic = PETSC_FALSE;
482b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
483b8840d6bSPeter Brune       if (i_r>qn->m-1) periodic = PETSC_TRUE;
4840c777b0cSPeter Brune     }
4850c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
4860c777b0cSPeter Brune     if (powell || periodic) {
4875ba6227bSPeter Brune       if (qn->monitor) {
4885ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
489516fe3c3SPeter 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);
4905ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4915ba6227bSPeter Brune       }
4925ba6227bSPeter Brune       i_r = -1;
4935ba6227bSPeter Brune       /* general purpose update */
4945ba6227bSPeter Brune       if (snes->ops->update) {
4955ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
4965ba6227bSPeter Brune       }
4970c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
4980ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
4998cba9625SJed Brown         ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
5000ecc9e1dSPeter Brune       }
5010ecc9e1dSPeter Brune     }
50270d3b23bSPeter Brune     /* general purpose update */
50370d3b23bSPeter Brune     if (snes->ops->update) {
50470d3b23bSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
50570d3b23bSPeter Brune     }
5065ba6227bSPeter Brune   }
5074b11644fSPeter Brune   if (i == snes->max_its) {
5084b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
5094b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
5104b11644fSPeter Brune   }
5114b11644fSPeter Brune   PetscFunctionReturn(0);
5124b11644fSPeter Brune }
5134b11644fSPeter Brune 
5144b11644fSPeter Brune #undef __FUNCT__
5154b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
5164b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
5174b11644fSPeter Brune {
5189f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
5194b11644fSPeter Brune   PetscErrorCode ierr;
520335fb975SPeter Brune 
5214b11644fSPeter Brune   PetscFunctionBegin;
522b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
523b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
5248e231d97SPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha,
5258e231d97SPeter Brune                       qn->m, PetscScalar, &qn->beta,
5268e231d97SPeter Brune                       qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr);
5278e231d97SPeter Brune 
52818aa0c0cSPeter Brune   if (qn->singlereduction) {
5298e231d97SPeter Brune     ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat,
5308e231d97SPeter Brune                         qn->m, PetscScalar, &qn->dFtdX,
5318e231d97SPeter Brune                         qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr);
5328e231d97SPeter Brune   }
533fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
534335fb975SPeter Brune 
535335fb975SPeter Brune   /* set up the line search */
5360c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5378e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
5388e231d97SPeter Brune   }
5394b11644fSPeter Brune   PetscFunctionReturn(0);
5404b11644fSPeter Brune }
5414b11644fSPeter Brune 
5424b11644fSPeter Brune #undef __FUNCT__
5434b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
5444b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
5454b11644fSPeter Brune {
5464b11644fSPeter Brune   PetscErrorCode ierr;
5479f83bee8SJed Brown   SNES_QN        *qn;
5480adebc6cSBarry Smith 
5494b11644fSPeter Brune   PetscFunctionBegin;
5504b11644fSPeter Brune   if (snes->data) {
5519f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
552b8840d6bSPeter Brune     if (qn->U) {
553b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
5544b11644fSPeter Brune     }
555b8840d6bSPeter Brune     if (qn->V) {
556b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
5574b11644fSPeter Brune     }
55818aa0c0cSPeter Brune     if (qn->singlereduction) {
5598e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
5608e231d97SPeter Brune     }
5618e231d97SPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr);
5624b11644fSPeter Brune   }
5634b11644fSPeter Brune   PetscFunctionReturn(0);
5644b11644fSPeter Brune }
5654b11644fSPeter Brune 
5664b11644fSPeter Brune #undef __FUNCT__
5674b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
5684b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
5694b11644fSPeter Brune {
5704b11644fSPeter Brune   PetscErrorCode ierr;
5716e111a19SKarl Rupp 
5724b11644fSPeter Brune   PetscFunctionBegin;
5734b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
5744b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
575bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
5764b11644fSPeter Brune   PetscFunctionReturn(0);
5774b11644fSPeter Brune }
5784b11644fSPeter Brune 
5794b11644fSPeter Brune #undef __FUNCT__
5804b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
5814b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
5824b11644fSPeter Brune {
5834b11644fSPeter Brune 
5844b11644fSPeter Brune   PetscErrorCode    ierr;
5852150357eSBarry Smith   SNES_QN           *qn    = (SNES_QN*)snes->data;
58688f769c5SPeter Brune   PetscBool         monflg = PETSC_FALSE,flg;
587f1c6b773SPeter Brune   SNESLineSearch    linesearch;
5882150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
5892150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
5902150357eSBarry Smith 
5914b11644fSPeter Brune   PetscFunctionBegin;
5924b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
5930298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
5940298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
5950298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr);
5960298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
5970298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
59888f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
59988f769c5SPeter Brune   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
60088f769c5SPeter Brune 
60188f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
60288f769c5SPeter Brune   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
60388f769c5SPeter Brune 
604b8840d6bSPeter Brune   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,
6050298fd71SBarry Smith                           (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr);
6064b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
6079e764e56SPeter Brune   if (!snes->linesearch) {
6087601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
6091a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
6109e764e56SPeter Brune   }
61144f7e39eSPeter Brune   if (monflg) {
612ce94432eSBarry Smith     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
61344f7e39eSPeter Brune   }
6144b11644fSPeter Brune   PetscFunctionReturn(0);
6154b11644fSPeter Brune }
6164b11644fSPeter Brune 
61792c02d66SPeter Brune 
6180c777b0cSPeter Brune #undef __FUNCT__
6190c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType"
6200c777b0cSPeter Brune /*@
6210c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
6220c777b0cSPeter Brune 
6230c777b0cSPeter Brune     Logically Collective on SNES
6240c777b0cSPeter Brune 
6250c777b0cSPeter Brune     Input Parameters:
6260c777b0cSPeter Brune +   snes - the iterative context
6270c777b0cSPeter Brune -   rtype - restart type
6280c777b0cSPeter Brune 
6290c777b0cSPeter Brune     Options Database:
6300c777b0cSPeter Brune +   -snes_qn_restart_type<powell,periodic,none> - set the restart type
631b8840d6bSPeter Brune -   -snes_qn_restart[10] - sets the number of iterations before restart for periodic
6320c777b0cSPeter Brune 
6330c777b0cSPeter Brune     Level: intermediate
6340c777b0cSPeter Brune 
6350c777b0cSPeter Brune     SNESQNRestartTypes:
6360c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
6370c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
6380c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
6390c777b0cSPeter Brune 
6400c777b0cSPeter Brune     Notes:
6410c777b0cSPeter Brune     The default line search used is the L2 line search and it requires two additional function evaluations.
6420c777b0cSPeter Brune 
6430c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch
6440c777b0cSPeter Brune @*/
6452150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
6462150357eSBarry Smith {
6470c777b0cSPeter Brune   PetscErrorCode ierr;
6486e111a19SKarl Rupp 
6490c777b0cSPeter Brune   PetscFunctionBegin;
6500c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6510c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
6520c777b0cSPeter Brune   PetscFunctionReturn(0);
6530c777b0cSPeter Brune }
6540c777b0cSPeter Brune 
6550c777b0cSPeter Brune #undef __FUNCT__
6560c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType"
6570c777b0cSPeter Brune /*@
6580c777b0cSPeter Brune     SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN.
6590c777b0cSPeter Brune 
6600c777b0cSPeter Brune     Logically Collective on SNES
6610c777b0cSPeter Brune 
6620c777b0cSPeter Brune     Input Parameters:
6630c777b0cSPeter Brune +   snes - the iterative context
6640c777b0cSPeter Brune -   stype - scale type
6650c777b0cSPeter Brune 
6660c777b0cSPeter Brune     Options Database:
6670c777b0cSPeter Brune .   -snes_qn_scale_type<shanno,none,linesearch,jacobian>
6680c777b0cSPeter Brune 
6690c777b0cSPeter Brune     Level: intermediate
6700c777b0cSPeter Brune 
6710c777b0cSPeter Brune     SNESQNSelectTypes:
6720c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
6730c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
6740c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
6750c777b0cSPeter Brune -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
6760c777b0cSPeter Brune 
6770c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
6780c777b0cSPeter Brune @*/
6790c777b0cSPeter Brune 
6802150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
6812150357eSBarry Smith {
6820c777b0cSPeter Brune   PetscErrorCode ierr;
6836e111a19SKarl Rupp 
6840c777b0cSPeter Brune   PetscFunctionBegin;
6850c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6860c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
6870c777b0cSPeter Brune   PetscFunctionReturn(0);
6880c777b0cSPeter Brune }
6890c777b0cSPeter Brune 
6900c777b0cSPeter Brune #undef __FUNCT__
6910c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN"
6920adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
6930adebc6cSBarry Smith {
6940c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
6956e111a19SKarl Rupp 
6960c777b0cSPeter Brune   PetscFunctionBegin;
6970c777b0cSPeter Brune   qn->scale_type = stype;
6980c777b0cSPeter Brune   PetscFunctionReturn(0);
6990c777b0cSPeter Brune }
7000c777b0cSPeter Brune 
7010c777b0cSPeter Brune #undef __FUNCT__
7020c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN"
7030adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
7040adebc6cSBarry Smith {
7050c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7066e111a19SKarl Rupp 
7070c777b0cSPeter Brune   PetscFunctionBegin;
7080c777b0cSPeter Brune   qn->restart_type = rtype;
7090c777b0cSPeter Brune   PetscFunctionReturn(0);
7100c777b0cSPeter Brune }
7110c777b0cSPeter Brune 
7124b11644fSPeter Brune /* -------------------------------------------------------------------------- */
7134b11644fSPeter Brune /*MC
7144b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
7154b11644fSPeter Brune 
7166cc8130cSPeter Brune       Options Database:
7176cc8130cSPeter Brune 
7186cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
7191867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
7201867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
72172835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
72244f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
7236cc8130cSPeter Brune 
724b8840d6bSPeter Brune       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
725b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
726b8840d6bSPeter Brune       updates.
7276cc8130cSPeter Brune 
7281867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
7291867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
7301867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
7311867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
7321867fe5bSPeter Brune 
7336cc8130cSPeter Brune       References:
7346cc8130cSPeter Brune 
7356cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
7366cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
7376cc8130cSPeter Brune 
738b8840d6bSPeter Brune       Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker
739b8840d6bSPeter Brune       SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
740b8840d6bSPeter Brune 
7414b11644fSPeter Brune 
7424b11644fSPeter Brune       Level: beginner
7434b11644fSPeter Brune 
74404d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
7456cc8130cSPeter Brune 
7464b11644fSPeter Brune M*/
7474b11644fSPeter Brune #undef __FUNCT__
7484b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
7498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
7504b11644fSPeter Brune {
7514b11644fSPeter Brune   PetscErrorCode ierr;
7529f83bee8SJed Brown   SNES_QN        *qn;
7534b11644fSPeter Brune 
7544b11644fSPeter Brune   PetscFunctionBegin;
7554b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
7564b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
7574b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
7584b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
7594b11644fSPeter Brune   snes->ops->view           = 0;
7604b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
7614b11644fSPeter Brune 
76247f26062SPeter Brune   snes->pcside = PC_LEFT;
76347f26062SPeter Brune   snes->functype = SNES_FUNCTION_PRECONDITIONED;
76447f26062SPeter Brune 
76542f4f86dSBarry Smith   snes->usespc  = PETSC_TRUE;
76642f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
76742f4f86dSBarry Smith 
76888976e71SPeter Brune   if (!snes->tolerancesset) {
7690e444f03SPeter Brune     snes->max_funcs = 30000;
7700e444f03SPeter Brune     snes->max_its   = 10000;
77188976e71SPeter Brune   }
7720e444f03SPeter Brune 
7739f83bee8SJed Brown   ierr                = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
7744b11644fSPeter Brune   snes->data          = (void*) qn;
7750ecc9e1dSPeter Brune   qn->m               = 10;
776b21d5a53SPeter Brune   qn->scaling         = 1.0;
7770298fd71SBarry Smith   qn->U               = NULL;
7780298fd71SBarry Smith   qn->V               = NULL;
7790298fd71SBarry Smith   qn->dXtdF           = NULL;
7800298fd71SBarry Smith   qn->dFtdX           = NULL;
7810298fd71SBarry Smith   qn->dXdFmat         = NULL;
7820298fd71SBarry Smith   qn->monitor         = NULL;
78318aa0c0cSPeter Brune   qn->singlereduction = PETSC_FALSE;
784b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
7850c777b0cSPeter Brune   qn->scale_type      = SNES_QN_SCALE_SHANNO;
7860c777b0cSPeter Brune   qn->restart_type    = SNES_QN_RESTART_POWELL;
787b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
788ea630c6eSPeter Brune 
789bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
790bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
7914b11644fSPeter Brune   PetscFunctionReturn(0);
7924b11644fSPeter Brune }
793