xref: /petsc/src/snes/impls/qn/qn.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h>
34b11644fSPeter Brune 
48e231d97SPeter Brune #define H(i,j)  qn->dXdFmat[i*qn->m + j]
58e231d97SPeter Brune 
61efc8c45SPeter Brune const char *const SNESQNScaleTypes[] =        {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
71efc8c45SPeter Brune const char *const SNESQNRestartTypes[] =      {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
860c08b40SPeter Brune const char *const SNESQNTypes[] =             {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0};
9b8840d6bSPeter Brune 
104b11644fSPeter Brune typedef struct {
11b8840d6bSPeter Brune   Vec               *U;                   /* Stored past states (vary from method to method) */
12b8840d6bSPeter Brune   Vec               *V;                   /* Stored past states (vary from method to method) */
138e231d97SPeter Brune   PetscInt          m;                    /* The number of kept previous steps */
145c7a0a03SPeter Brune   PetscReal         *lambda;              /* The line search history of the method */
155c7a0a03SPeter Brune   PetscReal         *norm;                /* norms of the steps */
168e231d97SPeter Brune   PetscScalar       *alpha, *beta;
178e231d97SPeter Brune   PetscScalar       *dXtdF, *dFtdX, *YtdX;
1818aa0c0cSPeter Brune   PetscBool         singlereduction;      /* Aggregated reduction implementation */
198e231d97SPeter Brune   PetscScalar       *dXdFmat;             /* A matrix of values for dX_i dot dF_j */
2044f7e39eSPeter Brune   PetscViewer       monitor;
216bf1b2e5SPeter Brune   PetscReal         powell_gamma;         /* Powell angle restart condition */
22b21d5a53SPeter Brune   PetscReal         scaling;              /* scaling of H0 */
23b8840d6bSPeter Brune   SNESQNType        type;                 /* the type of quasi-newton method used */
2488f769c5SPeter Brune   SNESQNScaleType   scale_type;           /* the type of scaling used */
250c777b0cSPeter Brune   SNESQNRestartType restart_type;         /* determine the frequency and type of restart conditions */
269f83bee8SJed Brown } SNES_QN;
274b11644fSPeter Brune 
285051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D)
292150357eSBarry Smith {
304b11644fSPeter Brune   PetscErrorCode     ierr;
319f83bee8SJed Brown   SNES_QN            *qn = (SNES_QN*)snes->data;
32b8840d6bSPeter Brune   Vec                W   = snes->work[3];
33b8840d6bSPeter Brune   Vec                *U  = qn->U;
34b8840d6bSPeter Brune   PetscInt           m = qn->m;
355c7a0a03SPeter Brune   PetscInt           k,i,j,lits,l = m;
365c7a0a03SPeter Brune   PetscReal          unorm,a,b;
375c7a0a03SPeter Brune   PetscReal          *lambda=qn->lambda;
38b8840d6bSPeter Brune   PetscScalar        gdot;
3959e7931aSPeter Brune   PetscReal          udot;
402150357eSBarry Smith 
41b8840d6bSPeter Brune   PetscFunctionBegin;
42b8840d6bSPeter Brune   if (it < m) l = it;
435c7a0a03SPeter Brune   if (it > 0) {
445c7a0a03SPeter Brune     k = (it-1)%l;
455c7a0a03SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr);
46a49fa0d8SPeter Brune     ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr);
475051edb1SPeter Brune     ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr);
485c7a0a03SPeter Brune     if (qn->monitor) {
495c7a0a03SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
506bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %D of %D by lambda: %14.12e \n",k,l,(double)lambda[k]);CHKERRQ(ierr);
515c7a0a03SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
525c7a0a03SPeter Brune     }
535c7a0a03SPeter Brune   }
54b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
55d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr);
56422a814eSBarry Smith     SNESCheckKSPSolve(snes);
57b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
58b8840d6bSPeter Brune     snes->linear_its += lits;
59b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
60b8840d6bSPeter Brune   } else {
61b8840d6bSPeter Brune     ierr = VecCopy(D,Y);CHKERRQ(ierr);
62b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
63b8840d6bSPeter Brune   }
64b8840d6bSPeter Brune 
65b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
66b8840d6bSPeter Brune   for (i = 0; i < l-1; i++) {
675c7a0a03SPeter Brune     j = (it+i-l)%l;
685c7a0a03SPeter Brune     k = (it+i-l+1)%l;
695c7a0a03SPeter Brune     ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr);
705c7a0a03SPeter Brune     ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr);
7159e7931aSPeter Brune     unorm *= unorm;
7259e7931aSPeter Brune     udot = PetscRealPart(gdot);
735c7a0a03SPeter Brune     a = (lambda[j]/lambda[k]);
745c7a0a03SPeter Brune     b = -(1.-lambda[j]);
7559e7931aSPeter Brune     a *= udot/unorm;
7659e7931aSPeter Brune     b *= udot/unorm;
775c7a0a03SPeter Brune     ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr);
785c7a0a03SPeter Brune 
79b8840d6bSPeter Brune     if (qn->monitor) {
80b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
816bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D and %D, gdot: %14.12e\n",k,j,(double)PetscRealPart(gdot));CHKERRQ(ierr);
82b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
83b8840d6bSPeter Brune     }
84b8840d6bSPeter Brune   }
85b8840d6bSPeter Brune   if (l > 0) {
86b8840d6bSPeter Brune     k = (it-1)%l;
87b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
885c7a0a03SPeter Brune     ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr);
8959e7931aSPeter Brune     unorm *= unorm;
9059e7931aSPeter Brune     udot = PetscRealPart(gdot);
9159e7931aSPeter Brune     a = unorm/(unorm-lambda[k]*udot);
9259e7931aSPeter Brune     b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
93b8840d6bSPeter Brune     if (qn->monitor) {
94b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
956bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D: a: %14.12e b: %14.12e \n",k,(double)a,(double)b);CHKERRQ(ierr);
96b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
97b8840d6bSPeter Brune     }
985c7a0a03SPeter Brune     ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr);
99b8840d6bSPeter Brune   }
1005c7a0a03SPeter Brune   l = m;
1015c7a0a03SPeter Brune   if (it+1<m)l=it+1;
1025c7a0a03SPeter Brune   k = it%l;
1035c7a0a03SPeter Brune   if (qn->monitor) {
1045c7a0a03SPeter Brune     ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1056bdcc836SBarry Smith     ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %D of %D\n",k,l);CHKERRQ(ierr);
1065c7a0a03SPeter Brune     ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1075c7a0a03SPeter Brune   }
108b8840d6bSPeter Brune   PetscFunctionReturn(0);
109b8840d6bSPeter Brune }
110b8840d6bSPeter Brune 
1110adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1120adebc6cSBarry Smith {
113b8840d6bSPeter Brune   PetscErrorCode ierr;
114b8840d6bSPeter Brune   SNES_QN        *qn = (SNES_QN*)snes->data;
115b8840d6bSPeter Brune   Vec            W   = snes->work[3];
116b8840d6bSPeter Brune   Vec            *U  = qn->U;
117b8840d6bSPeter Brune   Vec            *T  = qn->V;
118b8840d6bSPeter Brune 
11984c577daSBarry Smith   /* ksp thing for Jacobian scaling */
120c9e59058SPeter Brune   PetscInt           h,k,j,i,lits;
121b8840d6bSPeter Brune   PetscInt           m = qn->m;
1225051edb1SPeter Brune   PetscScalar        gdot,udot;
123b8840d6bSPeter Brune   PetscInt           l = m;
1240adebc6cSBarry Smith 
125b8840d6bSPeter Brune   PetscFunctionBegin;
126b8840d6bSPeter Brune   if (it < m) l = it;
127b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
128b8840d6bSPeter Brune   if (l > 0) {
129b8840d6bSPeter Brune     k    = (it-1)%l;
130c9e59058SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr);
131b8840d6bSPeter Brune     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
132b8840d6bSPeter Brune     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
1335051edb1SPeter Brune     ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr);
1345051edb1SPeter Brune     ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr);
135b8840d6bSPeter Brune   }
136b8840d6bSPeter Brune 
137b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
138d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
139422a814eSBarry Smith     SNESCheckKSPSolve(snes);
140b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
141b8840d6bSPeter Brune     snes->linear_its += lits;
142b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
143b8840d6bSPeter Brune   } else {
144b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
145b8840d6bSPeter Brune   }
1465051edb1SPeter Brune 
1475051edb1SPeter Brune   /* inward recursion starting at the first update and working forward */
1485051edb1SPeter Brune   if (l > 0) {
1495051edb1SPeter Brune     for (i = 0; i < l-1; i++) {
150c9e59058SPeter Brune       j    = (it+i-l)%l;
151c9e59058SPeter Brune       k    = (it+i-l+1)%l;
152c9e59058SPeter Brune       h    = (it-1)%l;
153c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr);
154c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr);
155c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr);
156c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr);
1575051edb1SPeter Brune       ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr);
158c9e59058SPeter Brune       ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr);
1595051edb1SPeter Brune       if (qn->monitor) {
1605051edb1SPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1616bdcc836SBarry Smith         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D gdot: %14.12e\n", it, k, (double)PetscRealPart(gdot));CHKERRQ(ierr);
1625051edb1SPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1635051edb1SPeter Brune       }
1645051edb1SPeter Brune     }
1655051edb1SPeter Brune   }
166b8840d6bSPeter Brune   PetscFunctionReturn(0);
167b8840d6bSPeter Brune }
168b8840d6bSPeter Brune 
1690adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1700adebc6cSBarry Smith {
171b8840d6bSPeter Brune   PetscErrorCode ierr;
172b8840d6bSPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
173b8840d6bSPeter Brune   Vec            W      = snes->work[3];
174b8840d6bSPeter Brune   Vec            *dX    = qn->U;
175b8840d6bSPeter Brune   Vec            *dF    = qn->V;
176bd052dfeSPeter Brune   PetscScalar    *alpha = qn->alpha;
177bd052dfeSPeter Brune   PetscScalar    *beta  = qn->beta;
1788e231d97SPeter Brune   PetscScalar    *dXtdF = qn->dXtdF;
179b8840d6bSPeter Brune   PetscScalar    *dFtdX = qn->dFtdX;
1808e231d97SPeter Brune   PetscScalar    *YtdX  = qn->YtdX;
181bd052dfeSPeter Brune 
18284c577daSBarry Smith   /* ksp thing for Jacobian scaling */
1838e231d97SPeter Brune   PetscInt           k,i,j,g,lits;
1844b11644fSPeter Brune   PetscInt           m = qn->m;
1854b11644fSPeter Brune   PetscScalar        t;
1864b11644fSPeter Brune   PetscInt           l = m;
1878e231d97SPeter Brune 
1884b11644fSPeter Brune   PetscFunctionBegin;
1895ba6227bSPeter Brune   if (it < m) l = it;
1901c6523dcSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
191b8840d6bSPeter Brune   if (it > 0) {
192b8840d6bSPeter Brune     k    = (it - 1) % l;
193b8840d6bSPeter Brune     ierr = VecCopy(D,dF[k]);CHKERRQ(ierr);
194b8840d6bSPeter Brune     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
195b8840d6bSPeter Brune     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
196b8840d6bSPeter Brune     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
19718aa0c0cSPeter Brune     if (qn->singlereduction) {
1981c6523dcSPeter Brune       ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
1991c6523dcSPeter Brune       ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
2001c6523dcSPeter Brune       ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr);
2011c6523dcSPeter Brune       ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
2021c6523dcSPeter Brune       ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
2031c6523dcSPeter Brune       ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr);
204b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
205b8840d6bSPeter Brune         H(k, j) = dFtdX[j];
206b8840d6bSPeter Brune         H(j, k) = dXtdF[j];
207b8840d6bSPeter Brune       }
208b8840d6bSPeter Brune       /* copy back over to make the computation of alpha and beta easier */
2091aa26658SKarl Rupp       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
210b8840d6bSPeter Brune     } else {
211b8840d6bSPeter Brune       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
21223d44fbcSPeter Brune     }
21323d44fbcSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
214b8840d6bSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
215b8840d6bSPeter Brune     }
216b8840d6bSPeter Brune   }
217b8840d6bSPeter Brune 
2184b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
2194b11644fSPeter Brune   for (i=0; i<l; i++) {
220b21d5a53SPeter Brune     k = (it-i-1)%l;
22118aa0c0cSPeter Brune     if (qn->singlereduction) {
2228e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
2238e231d97SPeter Brune       t = YtdX[k];
2248e231d97SPeter Brune       for (j=0; j<i; j++) {
2258e231d97SPeter Brune         g  = (it-j-1)%l;
226b3cfccb2SPeter Brune         t -= alpha[g]*H(k, g);
2278e231d97SPeter Brune       }
2288e231d97SPeter Brune       alpha[k] = t/H(k,k);
2298e231d97SPeter Brune     } else {
2303af51624SPeter Brune       ierr     = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
2318e231d97SPeter Brune       alpha[k] = t/dXtdF[k];
2328e231d97SPeter Brune     }
23344f7e39eSPeter Brune     if (qn->monitor) {
2343af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2356bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha:        %14.12e\n", it, k, (double)PetscRealPart(alpha[k]));CHKERRQ(ierr);
2363af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
23744f7e39eSPeter Brune     }
2386bf1b2e5SPeter Brune     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
2394b11644fSPeter Brune   }
2404b11644fSPeter Brune 
2410c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
242d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
243422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2440ecc9e1dSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2450ecc9e1dSPeter Brune     snes->linear_its += lits;
246b8840d6bSPeter Brune     ierr              = VecCopy(W, Y);CHKERRQ(ierr);
2470ecc9e1dSPeter Brune   } else {
248b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
2490ecc9e1dSPeter Brune   }
25018aa0c0cSPeter Brune   if (qn->singlereduction) {
251b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
2528e231d97SPeter Brune   }
2534b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
254bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
255b21d5a53SPeter Brune     k = (it + i - l) % l;
25618aa0c0cSPeter Brune     if (qn->singlereduction) {
2578e231d97SPeter Brune       t = YtdX[k];
2588e231d97SPeter Brune       for (j = 0; j < i; j++) {
2598e231d97SPeter Brune         g  = (it + j - l) % l;
260b3cfccb2SPeter Brune         t += (alpha[g] - beta[g])*H(g, k);
2618e231d97SPeter Brune       }
2628e231d97SPeter Brune       beta[k] = t / H(k, k);
2638e231d97SPeter Brune     } else {
2646bf1b2e5SPeter Brune       ierr    = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
2658e231d97SPeter Brune       beta[k] = t / dXtdF[k];
2668e231d97SPeter Brune     }
26722d28d08SBarry Smith     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr);
26844f7e39eSPeter Brune     if (qn->monitor) {
2693af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2706bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha - beta: %14.12e\n", it, k, (double)PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
2713af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
27244f7e39eSPeter Brune     }
2734b11644fSPeter Brune   }
2744b11644fSPeter Brune   PetscFunctionReturn(0);
2754b11644fSPeter Brune }
2764b11644fSPeter Brune 
2774b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
2784b11644fSPeter Brune {
2794b11644fSPeter Brune   PetscErrorCode       ierr;
2809f83bee8SJed Brown   SNES_QN              *qn = (SNES_QN*) snes->data;
28115f5eeeaSPeter Brune   Vec                  X,Xold;
2820a3905e1SPeter Brune   Vec                  F,W;
28347f26062SPeter Brune   Vec                  Y,D,Dold;
284b8840d6bSPeter Brune   PetscInt             i, i_r;
285335fb975SPeter Brune   PetscReal            fnorm,xnorm,ynorm,gnorm;
286422a814eSBarry Smith   SNESLineSearchReason lssucceed;
287422a814eSBarry Smith   PetscBool            powell,periodic;
2881c6523dcSPeter Brune   PetscScalar          DolddotD,DolddotDold;
289b7281c8aSPeter Brune   SNESConvergedReason  reason;
2900ecc9e1dSPeter Brune 
29184c577daSBarry Smith   /* basically just a regular newton's method except for the application of the Jacobian */
2924b11644fSPeter Brune 
2936e111a19SKarl Rupp   PetscFunctionBegin;
2946c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
295c579b300SPatrick Farrell 
296fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
2979f3a0142SPeter Brune   F    = snes->vec_func;                /* residual vector */
2983af51624SPeter Brune   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
2990a3905e1SPeter Brune   W    = snes->work[3];
300b8840d6bSPeter Brune   X    = snes->vec_sol;                 /* solution vector */
301335fb975SPeter Brune   Xold = snes->work[0];
3023af51624SPeter Brune 
3033af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
304335fb975SPeter Brune   D    = snes->work[1];
305335fb975SPeter Brune   Dold = snes->work[2];
3064b11644fSPeter Brune 
3074b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
3084b11644fSPeter Brune 
309e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3104b11644fSPeter Brune   snes->iter = 0;
3114b11644fSPeter Brune   snes->norm = 0.;
312e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
31347f26062SPeter Brune 
314*efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
315be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
316*efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
317b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
318b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
319b7281c8aSPeter Brune       PetscFunctionReturn(0);
320b7281c8aSPeter Brune     }
32147f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
32247f26062SPeter Brune   } else {
323e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
32415f5eeeaSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3251aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
326c1c75074SPeter Brune 
32747f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
328422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
32947f26062SPeter Brune   }
330*efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
331be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
332*efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&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   } else {
33847f26062SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
33947f26062SPeter Brune   }
340b28a06ddSPeter Brune 
341e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3424b11644fSPeter Brune   snes->norm = fnorm;
343e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
344a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
3454b11644fSPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
3464b11644fSPeter Brune 
3474b11644fSPeter Brune   /* test convergence */
3484b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3494b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
35070d3b23bSPeter Brune 
351*efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_RIGHT) {
352*efd4aadfSBarry Smith     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
353*efd4aadfSBarry Smith     ierr = SNESSolve(snes->npc,snes->vec_rhs,X);CHKERRQ(ierr);
354*efd4aadfSBarry Smith     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
355*efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
356ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
357ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
358ddd40ce5SPeter Brune       PetscFunctionReturn(0);
359ddd40ce5SPeter Brune     }
360be95d8f1SBarry Smith     ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
3613cf07b75SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
3623cf07b75SPeter Brune   }
3633cf07b75SPeter Brune 
364f8e15203SPeter Brune   /* scale the initial update */
3650c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
366d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
367aeb49b86SAsbjørn Nilsen Riseth     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
3680ecc9e1dSPeter Brune   }
3690ecc9e1dSPeter Brune 
3705ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
3710a3905e1SPeter Brune     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
3720a3905e1SPeter Brune       PetscScalar ff,xf;
3730a3905e1SPeter Brune       ierr = VecCopy(Dold,Y);CHKERRQ(ierr);
3740a3905e1SPeter Brune       ierr = VecCopy(Xold,W);CHKERRQ(ierr);
3750a3905e1SPeter Brune       ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr);
3760a3905e1SPeter Brune       ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
3770a3905e1SPeter Brune       ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr);
3780a3905e1SPeter Brune       ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr);
3790a3905e1SPeter Brune       ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr);
3800a3905e1SPeter Brune       ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr);
3810a3905e1SPeter Brune       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
3826bdcc836SBarry Smith       if (qn->monitor) {
3836bdcc836SBarry Smith         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3846bdcc836SBarry Smith         ierr = PetscViewerASCIIPrintf(qn->monitor, "Shanno scaling %D %g\n", i,(double)qn->scaling);CHKERRQ(ierr);
3856bdcc836SBarry Smith         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3866bdcc836SBarry Smith       }
3870a3905e1SPeter Brune     }
388b8840d6bSPeter Brune     switch (qn->type) {
389b8840d6bSPeter Brune     case SNES_QN_BADBROYDEN:
390b8840d6bSPeter Brune       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
391b8840d6bSPeter Brune       break;
392b8840d6bSPeter Brune     case SNES_QN_BROYDEN:
3935051edb1SPeter Brune       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr);
394b8840d6bSPeter Brune       break;
395b8840d6bSPeter Brune     case SNES_QN_LBFGS:
396b8840d6bSPeter Brune       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
397b8840d6bSPeter Brune       break;
398b8840d6bSPeter Brune     }
39970d3b23bSPeter Brune     /* line search for lambda */
40070d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
4010a3905e1SPeter Brune     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
40215f5eeeaSPeter Brune     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
403f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
4049f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
405422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr);
406422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
407422a814eSBarry Smith     if (lssucceed) {
4089f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
4099f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
4109f3a0142SPeter Brune         break;
4119f3a0142SPeter Brune       }
4129f3a0142SPeter Brune     }
4130c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
41405ee83feSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
4150ecc9e1dSPeter Brune     }
4163af51624SPeter Brune 
41788d374b2SPeter Brune     /* convergence monitoring */
418b21d5a53SPeter 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);
419b21d5a53SPeter Brune 
420*efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_RIGHT) {
421*efd4aadfSBarry Smith       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
422*efd4aadfSBarry Smith       ierr = SNESSolve(snes->npc,snes->vec_rhs,X);CHKERRQ(ierr);
423*efd4aadfSBarry Smith       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
424*efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
425ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
426ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
427ddd40ce5SPeter Brune         PetscFunctionReturn(0);
428ddd40ce5SPeter Brune       }
429be95d8f1SBarry Smith       ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
430b28a06ddSPeter Brune     }
431b28a06ddSPeter Brune 
432360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
43371dbe336SPeter Brune     snes->norm = fnorm;
434360c497dSPeter Brune 
435a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
4368409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
4374b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
438d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4394b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
440*efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
441be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
442*efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
443b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
444b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
445b7281c8aSPeter Brune         PetscFunctionReturn(0);
446b7281c8aSPeter Brune       }
44788d374b2SPeter Brune     } else {
44888d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
44988d374b2SPeter Brune     }
4500c777b0cSPeter Brune     powell = PETSC_FALSE;
4516bdcc836SBarry Smith     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
4526bdcc836SBarry Smith       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
45323c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
45423c918e7SPeter Brune         ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr);
45523c918e7SPeter Brune       } else {
45623c918e7SPeter Brune         ierr = VecCopy(Dold,W);CHKERRQ(ierr);
45723c918e7SPeter Brune       }
45823c918e7SPeter Brune       ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr);
45923c918e7SPeter Brune       ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr);
46023c918e7SPeter Brune       ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr);
46123c918e7SPeter Brune       ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr);
4620c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
4630c777b0cSPeter Brune     }
4640c777b0cSPeter Brune     periodic = PETSC_FALSE;
465b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
466b8840d6bSPeter Brune       if (i_r>qn->m-1) periodic = PETSC_TRUE;
4670c777b0cSPeter Brune     }
4680c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
4690c777b0cSPeter Brune     if (powell || periodic) {
4705ba6227bSPeter Brune       if (qn->monitor) {
4715ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4726bdcc836SBarry Smith         if (powell) {
4736bdcc836SBarry Smith           ierr = PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %D\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold),i_r);CHKERRQ(ierr);
4746bdcc836SBarry Smith         } else {
4756bdcc836SBarry Smith           ierr = PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %D\n", i_r);CHKERRQ(ierr);
4766bdcc836SBarry Smith         }
4775ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4785ba6227bSPeter Brune       }
4795ba6227bSPeter Brune       i_r = -1;
4805ba6227bSPeter Brune       /* general purpose update */
4815ba6227bSPeter Brune       if (snes->ops->update) {
4825ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
4835ba6227bSPeter Brune       }
4840c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
485d1e9a80fSBarry Smith         ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
4860ecc9e1dSPeter Brune       }
4870ecc9e1dSPeter Brune     }
48870d3b23bSPeter Brune     /* general purpose update */
48970d3b23bSPeter Brune     if (snes->ops->update) {
49070d3b23bSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
49170d3b23bSPeter Brune     }
4925ba6227bSPeter Brune   }
4934b11644fSPeter Brune   if (i == snes->max_its) {
4944b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
4954b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
4964b11644fSPeter Brune   }
4974b11644fSPeter Brune   PetscFunctionReturn(0);
4984b11644fSPeter Brune }
4994b11644fSPeter Brune 
5004b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
5014b11644fSPeter Brune {
5029f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
5034b11644fSPeter Brune   PetscErrorCode ierr;
504fced5a79SAsbjørn Nilsen Riseth   DM             dm;
505335fb975SPeter Brune 
5064b11644fSPeter Brune   PetscFunctionBegin;
507fced5a79SAsbjørn Nilsen Riseth 
508fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
509fced5a79SAsbjørn Nilsen Riseth     ierr             = SNESGetDM(snes,&dm);CHKERRQ(ierr);
510fced5a79SAsbjørn Nilsen Riseth     ierr             = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr);
511fced5a79SAsbjørn Nilsen Riseth   }
512fced5a79SAsbjørn Nilsen Riseth 
513b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
51459e7931aSPeter Brune   if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
515dcca6d9dSJed Brown   ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr);
5168e231d97SPeter Brune 
51718aa0c0cSPeter Brune   if (qn->singlereduction) {
518dcca6d9dSJed Brown     ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr);
5198e231d97SPeter Brune   }
520fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
52160c08b40SPeter Brune   /* set method defaults */
5221efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
52360c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
52460c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
52560c08b40SPeter Brune     } else {
52660c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_SHANNO;
52760c08b40SPeter Brune     }
52860c08b40SPeter Brune   }
5291efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
53060c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
53160c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
53260c08b40SPeter Brune     } else {
53360c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
53460c08b40SPeter Brune     }
53560c08b40SPeter Brune   }
53660c08b40SPeter Brune 
5370c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5388e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
5398e231d97SPeter Brune   }
540*efd4aadfSBarry Smith   if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
5414b11644fSPeter Brune   PetscFunctionReturn(0);
5424b11644fSPeter Brune }
5434b11644fSPeter Brune 
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     }
5615c7a0a03SPeter Brune     ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr);
5624b11644fSPeter Brune   }
5634b11644fSPeter Brune   PetscFunctionReturn(0);
5644b11644fSPeter Brune }
5654b11644fSPeter Brune 
5664b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
5674b11644fSPeter Brune {
5684b11644fSPeter Brune   PetscErrorCode ierr;
5696e111a19SKarl Rupp 
5704b11644fSPeter Brune   PetscFunctionBegin;
5714b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
5724b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
573bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
5744b11644fSPeter Brune   PetscFunctionReturn(0);
5754b11644fSPeter Brune }
5764b11644fSPeter Brune 
5774416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(PetscOptionItems *PetscOptionsObject,SNES snes)
5784b11644fSPeter Brune {
5794b11644fSPeter Brune 
5804b11644fSPeter Brune   PetscErrorCode    ierr;
5812150357eSBarry Smith   SNES_QN           *qn    = (SNES_QN*)snes->data;
58288f769c5SPeter Brune   PetscBool         monflg = PETSC_FALSE,flg;
583f1c6b773SPeter Brune   SNESLineSearch    linesearch;
5842150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
5852150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
5861efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
5872150357eSBarry Smith 
5884b11644fSPeter Brune   PetscFunctionBegin;
589e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES QN options");CHKERRQ(ierr);
5900298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
5910298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
5920298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
5930298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
59488f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
59588f769c5SPeter Brune   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
59688f769c5SPeter Brune 
59788f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
59888f769c5SPeter Brune   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
59988f769c5SPeter Brune 
6000fdccdaeSBarry Smith   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr);
6011efc8c45SPeter Brune   if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);}
6024b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
6039e764e56SPeter Brune   if (!snes->linesearch) {
6047601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
60560c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
6061a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
60760c08b40SPeter Brune     } else if (qn->type == SNES_QN_BROYDEN) {
60860c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr);
60960c08b40SPeter Brune     } else {
61060c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
61160c08b40SPeter Brune     }
6129e764e56SPeter Brune   }
61344f7e39eSPeter Brune   if (monflg) {
614ce94432eSBarry Smith     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
61544f7e39eSPeter Brune   }
6164b11644fSPeter Brune   PetscFunctionReturn(0);
6174b11644fSPeter Brune }
6184b11644fSPeter Brune 
6195cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
6205cd83419SPeter Brune {
6215cd83419SPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
6225cd83419SPeter Brune   PetscBool      iascii;
6235cd83419SPeter Brune   PetscErrorCode ierr;
6245cd83419SPeter Brune 
6255cd83419SPeter Brune   PetscFunctionBegin;
6265cd83419SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
6275cd83419SPeter Brune   if (iascii) {
628*efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, restart type is %s, scale type is %s\n",SNESQNTypes[qn->type],SNESQNRestartTypes[qn->restart_type],SNESQNScaleTypes[qn->scale_type]);CHKERRQ(ierr);
6296bdcc836SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %D\n", qn->m);CHKERRQ(ierr);
6305cd83419SPeter Brune     if (qn->singlereduction) {
6315cd83419SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");CHKERRQ(ierr);
6325cd83419SPeter Brune     }
6335cd83419SPeter Brune   }
6345cd83419SPeter Brune   PetscFunctionReturn(0);
6355cd83419SPeter Brune }
63692c02d66SPeter Brune 
6370c777b0cSPeter Brune /*@
6380c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
6390c777b0cSPeter Brune 
6400c777b0cSPeter Brune     Logically Collective on SNES
6410c777b0cSPeter Brune 
6420c777b0cSPeter Brune     Input Parameters:
6430c777b0cSPeter Brune +   snes - the iterative context
6440c777b0cSPeter Brune -   rtype - restart type
6450c777b0cSPeter Brune 
6460c777b0cSPeter Brune     Options Database:
6470c777b0cSPeter Brune +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
64884c577daSBarry Smith -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
6490c777b0cSPeter Brune 
6500c777b0cSPeter Brune     Level: intermediate
6510c777b0cSPeter Brune 
6520c777b0cSPeter Brune     SNESQNRestartTypes:
6530c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
6540c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
6550c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
6560c777b0cSPeter Brune 
65784c577daSBarry Smith .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType
6580c777b0cSPeter Brune @*/
6592150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
6602150357eSBarry Smith {
6610c777b0cSPeter Brune   PetscErrorCode ierr;
6626e111a19SKarl Rupp 
6630c777b0cSPeter Brune   PetscFunctionBegin;
6640c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6650c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
6660c777b0cSPeter Brune   PetscFunctionReturn(0);
6670c777b0cSPeter Brune }
6680c777b0cSPeter Brune 
6690c777b0cSPeter Brune /*@
67084c577daSBarry Smith     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.
6710c777b0cSPeter Brune 
6720c777b0cSPeter Brune     Logically Collective on SNES
6730c777b0cSPeter Brune 
6740c777b0cSPeter Brune     Input Parameters:
6750c777b0cSPeter Brune +   snes - the iterative context
6760c777b0cSPeter Brune -   stype - scale type
6770c777b0cSPeter Brune 
6780c777b0cSPeter Brune     Options Database:
6790c777b0cSPeter Brune .   -snes_qn_scale_type <shanno,none,linesearch,jacobian>
6800c777b0cSPeter Brune 
6810c777b0cSPeter Brune     Level: intermediate
6820c777b0cSPeter Brune 
68384c577daSBarry Smith     SNESQNScaleTypes:
6840c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
6850c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
6860c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
687a01a0525SBarry Smith -   SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
688a01a0525SBarry Smith                              of QN and at ever restart.
6890c777b0cSPeter Brune 
690a01a0525SBarry Smith .keywords: scaling type
691a01a0525SBarry Smith 
692a01a0525SBarry Smith .seealso: SNES, SNESQN, SNESLineSearch, SNESQNScaleType, SNESetJacobian()
6930c777b0cSPeter Brune @*/
6940c777b0cSPeter Brune 
6952150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
6962150357eSBarry Smith {
6970c777b0cSPeter Brune   PetscErrorCode ierr;
6986e111a19SKarl Rupp 
6990c777b0cSPeter Brune   PetscFunctionBegin;
7000c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7010c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
7020c777b0cSPeter Brune   PetscFunctionReturn(0);
7030c777b0cSPeter Brune }
7040c777b0cSPeter Brune 
7050adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
7060adebc6cSBarry Smith {
7070c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7086e111a19SKarl Rupp 
7090c777b0cSPeter Brune   PetscFunctionBegin;
7100c777b0cSPeter Brune   qn->scale_type = stype;
7110c777b0cSPeter Brune   PetscFunctionReturn(0);
7120c777b0cSPeter Brune }
7130c777b0cSPeter Brune 
7140adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
7150adebc6cSBarry Smith {
7160c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7176e111a19SKarl Rupp 
7180c777b0cSPeter Brune   PetscFunctionBegin;
7190c777b0cSPeter Brune   qn->restart_type = rtype;
7200c777b0cSPeter Brune   PetscFunctionReturn(0);
7210c777b0cSPeter Brune }
7220c777b0cSPeter Brune 
7231efc8c45SPeter Brune /*@
7241efc8c45SPeter Brune     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
7251efc8c45SPeter Brune 
7261efc8c45SPeter Brune     Logically Collective on SNES
7271efc8c45SPeter Brune 
7281efc8c45SPeter Brune     Input Parameters:
7291efc8c45SPeter Brune +   snes - the iterative context
7301efc8c45SPeter Brune -   qtype - variant type
7311efc8c45SPeter Brune 
7321efc8c45SPeter Brune     Options Database:
73384c577daSBarry Smith .   -snes_qn_type <lbfgs,broyden,badbroyden>
7341efc8c45SPeter Brune 
7351efc8c45SPeter Brune     Level: beginner
7361efc8c45SPeter Brune 
7371efc8c45SPeter Brune     SNESQNTypes:
7381efc8c45SPeter Brune +   SNES_QN_LBFGS - LBFGS variant
7391efc8c45SPeter Brune .   SNES_QN_BROYDEN - Broyden variant
7401efc8c45SPeter Brune -   SNES_QN_BADBROYDEN - Bad Broyden variant
7411efc8c45SPeter Brune 
74284c577daSBarry Smith .keywords: SNES, SNESQN, type, set, SNESQNType
7431efc8c45SPeter Brune @*/
7441efc8c45SPeter Brune 
7451efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
7461efc8c45SPeter Brune {
7471efc8c45SPeter Brune   PetscErrorCode ierr;
7481efc8c45SPeter Brune 
7491efc8c45SPeter Brune   PetscFunctionBegin;
7501efc8c45SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7511efc8c45SPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr);
7521efc8c45SPeter Brune   PetscFunctionReturn(0);
7531efc8c45SPeter Brune }
7541efc8c45SPeter Brune 
7551efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
7561efc8c45SPeter Brune {
7571efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7581efc8c45SPeter Brune 
7591efc8c45SPeter Brune   PetscFunctionBegin;
7601efc8c45SPeter Brune   qn->type = qtype;
7611efc8c45SPeter Brune   PetscFunctionReturn(0);
7621efc8c45SPeter Brune }
7631efc8c45SPeter Brune 
7644b11644fSPeter Brune /* -------------------------------------------------------------------------- */
7654b11644fSPeter Brune /*MC
7664b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
7674b11644fSPeter Brune 
7686cc8130cSPeter Brune       Options Database:
7696cc8130cSPeter Brune 
77084c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
77184c577daSBarry Smith +     -snes_qn_restart_type <powell,periodic,none> - set the restart type
7726bdcc836SBarry Smith .     -snes_qn_powell_gamma - Angle condition for restart.
7731867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
77484c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
77584c577daSBarry Smith .     -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian
77672835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
77784c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
7786cc8130cSPeter Brune 
779b8840d6bSPeter Brune       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
780b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
781b8840d6bSPeter Brune       updates.
7826cc8130cSPeter Brune 
7831867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
7841867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
78584c577daSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
78684c577daSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
7871867fe5bSPeter Brune 
7882d547940SBarry Smith       Uses left nonlinear preconditioning by default.
7892d547940SBarry Smith 
7906cc8130cSPeter Brune       References:
79196a0c994SBarry Smith +   1. -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
79296a0c994SBarry Smith .   2. -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
7930a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
79496a0c994SBarry Smith .   3. -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
7950a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
79696a0c994SBarry Smith -   4. -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
7974f02bc6aSBarry Smith        SIAM Review, 57(4), 2015
7984b11644fSPeter Brune 
7994b11644fSPeter Brune       Level: beginner
8004b11644fSPeter Brune 
80104d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
8026cc8130cSPeter Brune 
8034b11644fSPeter Brune M*/
8048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
8054b11644fSPeter Brune {
8064b11644fSPeter Brune   PetscErrorCode ierr;
8079f83bee8SJed Brown   SNES_QN        *qn;
8084b11644fSPeter Brune 
8094b11644fSPeter Brune   PetscFunctionBegin;
8104b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
8114b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
8124b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
8134b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
8145cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
8154b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
8164b11644fSPeter Brune 
817*efd4aadfSBarry Smith   snes->npcside= PC_LEFT;
81847f26062SPeter Brune 
819*efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
82042f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
82142f4f86dSBarry Smith 
8224fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
8234fc747eaSLawrence Mitchell 
82488976e71SPeter Brune   if (!snes->tolerancesset) {
8250e444f03SPeter Brune     snes->max_funcs = 30000;
8260e444f03SPeter Brune     snes->max_its   = 10000;
82788976e71SPeter Brune   }
8280e444f03SPeter Brune 
829b00a9115SJed Brown   ierr                = PetscNewLog(snes,&qn);CHKERRQ(ierr);
8304b11644fSPeter Brune   snes->data          = (void*) qn;
8310ecc9e1dSPeter Brune   qn->m               = 10;
832b21d5a53SPeter Brune   qn->scaling         = 1.0;
8330298fd71SBarry Smith   qn->U               = NULL;
8340298fd71SBarry Smith   qn->V               = NULL;
8350298fd71SBarry Smith   qn->dXtdF           = NULL;
8360298fd71SBarry Smith   qn->dFtdX           = NULL;
8370298fd71SBarry Smith   qn->dXdFmat         = NULL;
8380298fd71SBarry Smith   qn->monitor         = NULL;
8391c6523dcSPeter Brune   qn->singlereduction = PETSC_TRUE;
840b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
84160c08b40SPeter Brune   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
84260c08b40SPeter Brune   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
843b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
844ea630c6eSPeter Brune 
845bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
846bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
8471efc8c45SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr);
8484b11644fSPeter Brune   PetscFunctionReturn(0);
8494b11644fSPeter Brune }
850