xref: /petsc/src/snes/impls/qn/qn.c (revision 4b11644f5537bbf1017da711cdde6651531d3f29)
1*4b11644fSPeter Brune 
2*4b11644fSPeter Brune #include <private/snesimpl.h>
3*4b11644fSPeter Brune 
4*4b11644fSPeter Brune typedef struct {
5*4b11644fSPeter Brune   PetscScalar lambda;   /* The default step length for the update */
6*4b11644fSPeter Brune   Vec * dX;  /* The change in X */
7*4b11644fSPeter Brune   Vec * dF;  /* The change in F */
8*4b11644fSPeter Brune   PetscInt m;     /* the number of kept previous steps */
9*4b11644fSPeter Brune } QNContext;
10*4b11644fSPeter Brune 
11*4b11644fSPeter Brune #undef __FUNCT__
12*4b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private"
13*4b11644fSPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) {
14*4b11644fSPeter Brune 
15*4b11644fSPeter Brune   PetscErrorCode ierr;
16*4b11644fSPeter Brune 
17*4b11644fSPeter Brune   QNContext * qn = (QNContext *)snes->data;
18*4b11644fSPeter Brune 
19*4b11644fSPeter Brune   Vec * dX = qn->dX;
20*4b11644fSPeter Brune   Vec * dF = qn->dF;
21*4b11644fSPeter Brune 
22*4b11644fSPeter Brune   PetscInt k, i;
23*4b11644fSPeter Brune   PetscInt m = qn->m;
24*4b11644fSPeter Brune   PetscScalar * alpha = PETSC_NULL;
25*4b11644fSPeter Brune   PetscScalar * beta = PETSC_NULL;
26*4b11644fSPeter Brune   PetscScalar * rho = PETSC_NULL;
27*4b11644fSPeter Brune   PetscScalar t;
28*4b11644fSPeter Brune   PetscInt l = m;
29*4b11644fSPeter Brune 
30*4b11644fSPeter Brune   PetscFunctionBegin;
31*4b11644fSPeter Brune 
32*4b11644fSPeter Brune   if (it < m) l = it+1;
33*4b11644fSPeter Brune   ierr = PetscMalloc3(m, PetscScalar, &alpha, m, PetscScalar, &beta, m, PetscScalar, &rho);CHKERRQ(ierr);
34*4b11644fSPeter Brune 
35*4b11644fSPeter Brune   /* precalculate alpha, beta, rho corresponding to the normal indices*/
36*4b11644fSPeter Brune   for (i = 0; i < l; i++) {
37*4b11644fSPeter Brune     ierr = VecDot(dX[i], dF[i], &t);CHKERRQ(ierr);
38*4b11644fSPeter Brune     rho[i] = 1. / t;
39*4b11644fSPeter Brune     beta[i] = 0;
40*4b11644fSPeter Brune     alpha[i] = 0;
41*4b11644fSPeter Brune   }
42*4b11644fSPeter Brune 
43*4b11644fSPeter Brune   ierr = VecCopy(g, z);CHKERRQ(ierr);
44*4b11644fSPeter Brune 
45*4b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
46*4b11644fSPeter Brune   for (i = 0; i < l; i++) {
47*4b11644fSPeter Brune     k = (it - i) % m;
48*4b11644fSPeter Brune     ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr);
49*4b11644fSPeter Brune     alpha[k] = t / rho[k];
50*4b11644fSPeter Brune     ierr = PetscPrintf(PETSC_COMM_WORLD, " %d: %e ", k, -alpha[k]);CHKERRQ(ierr);
51*4b11644fSPeter Brune     ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr);
52*4b11644fSPeter Brune   }
53*4b11644fSPeter Brune 
54*4b11644fSPeter Brune   /* inner application of the initial inverse jacobian approximation */
55*4b11644fSPeter Brune   /* right now it's just the identity. Nothing needs to go here. */
56*4b11644fSPeter Brune 
57*4b11644fSPeter Brune   /* inward recursion starting at the first update and working forward*/
58*4b11644fSPeter Brune   for (i = l - 1; i >= 0; i--) {
59*4b11644fSPeter Brune     k = (it - i) % m;
60*4b11644fSPeter Brune     ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr);
61*4b11644fSPeter Brune     beta[k] = rho[k]*t;
62*4b11644fSPeter Brune     ierr = PetscPrintf(PETSC_COMM_WORLD, " %d: %e %e", k, alpha[k] - beta[k], rho[k]);CHKERRQ(ierr);
63*4b11644fSPeter Brune     ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]);
64*4b11644fSPeter Brune   }
65*4b11644fSPeter Brune   ierr = PetscPrintf(PETSC_COMM_WORLD, "\n", k);CHKERRQ(ierr);
66*4b11644fSPeter Brune   ierr = PetscFree3(alpha, beta, rho);CHKERRQ(ierr);
67*4b11644fSPeter Brune 
68*4b11644fSPeter Brune   PetscFunctionReturn(0);
69*4b11644fSPeter Brune }
70*4b11644fSPeter Brune 
71*4b11644fSPeter Brune #undef __FUNCT__
72*4b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
73*4b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
74*4b11644fSPeter Brune {
75*4b11644fSPeter Brune 
76*4b11644fSPeter Brune   PetscErrorCode ierr;
77*4b11644fSPeter Brune   QNContext * qn = (QNContext*) snes->data;
78*4b11644fSPeter Brune 
79*4b11644fSPeter Brune   Vec x;
80*4b11644fSPeter Brune   Vec f;
81*4b11644fSPeter Brune   Vec p, pold;
82*4b11644fSPeter Brune 
83*4b11644fSPeter Brune   PetscInt i, j, k, l;
84*4b11644fSPeter Brune 
85*4b11644fSPeter Brune   PetscReal fnorm;
86*4b11644fSPeter Brune   PetscScalar gdot;
87*4b11644fSPeter Brune   PetscInt m = qn->m;
88*4b11644fSPeter Brune 
89*4b11644fSPeter Brune   Vec * W = qn->dX;
90*4b11644fSPeter Brune   Vec * V = qn->dF;
91*4b11644fSPeter Brune 
92*4b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
93*4b11644fSPeter Brune   PetscFunctionBegin;
94*4b11644fSPeter Brune 
95*4b11644fSPeter Brune   x = snes->vec_sol;
96*4b11644fSPeter Brune   f = snes->vec_func;
97*4b11644fSPeter Brune   p = snes->vec_sol_update;
98*4b11644fSPeter Brune   pold = snes->work[0];
99*4b11644fSPeter Brune 
100*4b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
101*4b11644fSPeter Brune 
102*4b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
103*4b11644fSPeter Brune   snes->iter = 0;
104*4b11644fSPeter Brune   snes->norm = 0.;
105*4b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
106*4b11644fSPeter Brune   ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
107*4b11644fSPeter Brune   if (snes->domainerror) {
108*4b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
109*4b11644fSPeter Brune     PetscFunctionReturn(0);
110*4b11644fSPeter Brune   }
111*4b11644fSPeter Brune   ierr = VecNorm(f, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
112*4b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
113*4b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
114*4b11644fSPeter Brune   snes->norm = fnorm;
115*4b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
116*4b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
117*4b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
118*4b11644fSPeter Brune 
119*4b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
120*4b11644fSPeter Brune   snes->ttol = fnorm*snes->rtol;
121*4b11644fSPeter Brune   /* test convergence */
122*4b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
123*4b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
124*4b11644fSPeter Brune   ierr = VecCopy(f, pold);CHKERRQ(ierr);
125*4b11644fSPeter Brune   ierr = VecAXPY(x, -1.0, pold);CHKERRQ(ierr);
126*4b11644fSPeter Brune   for(i = 0; i < snes->max_its; i++) {
127*4b11644fSPeter Brune 
128*4b11644fSPeter Brune     ierr = SNESComputeFunction(snes, x, f);CHKERRQ(ierr);
129*4b11644fSPeter Brune     ierr = VecNorm(f, NORM_2, &fnorm);CHKERRQ(ierr);
130*4b11644fSPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
131*4b11644fSPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
132*4b11644fSPeter Brune     snes->norm = fnorm;
133*4b11644fSPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
134*4b11644fSPeter Brune     SNESLogConvHistory(snes,fnorm,i+1);
135*4b11644fSPeter Brune     ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
136*4b11644fSPeter Brune 
137*4b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
138*4b11644fSPeter Brune     snes->ttol = fnorm*snes->rtol;
139*4b11644fSPeter Brune     /* test convergence */
140*4b11644fSPeter Brune     ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
141*4b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
142*4b11644fSPeter Brune     k = (i) % m;
143*4b11644fSPeter Brune     l = m;
144*4b11644fSPeter Brune     if (i < l) l = i;
145*4b11644fSPeter Brune     ierr = VecCopy(f, p);CHKERRQ(ierr);
146*4b11644fSPeter Brune     for (j=0; j<k; j++) {                                     /* p = product_{j<i} [I+v(j)w(j)^T]*p */
147*4b11644fSPeter Brune       ierr = VecDot(W[j],p,&gdot);CHKERRQ(ierr);
148*4b11644fSPeter Brune       ierr = VecAXPY(p,gdot,V[j]);CHKERRQ(ierr);
149*4b11644fSPeter Brune     }
150*4b11644fSPeter Brune     ierr = VecCopy(pold,W[k]);CHKERRQ(ierr);                  /* w[i] = pold   */
151*4b11644fSPeter Brune     ierr = VecAXPY(pold,-1.0,p);CHKERRQ(ierr);                /* v[i] =         p         */
152*4b11644fSPeter Brune     ierr = VecDot(W[k],pold,&gdot);CHKERRQ(ierr);             /*        ----------------- */
153*4b11644fSPeter Brune     ierr = VecCopy(p,V[k]);CHKERRQ(ierr);                     /*         w[i]'*(Pold - p) */
154*4b11644fSPeter Brune     ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr);
155*4b11644fSPeter Brune 
156*4b11644fSPeter Brune     ierr = VecDot(W[k],p,&gdot);CHKERRQ(ierr);                /* p = (I + v[i]*w[i]')*p   */
157*4b11644fSPeter Brune     ierr = VecAXPY(p,gdot,V[k]);CHKERRQ(ierr);
158*4b11644fSPeter Brune     ierr = VecCopy(p,pold);CHKERRQ(ierr);
159*4b11644fSPeter Brune 
160*4b11644fSPeter Brune     ierr = PetscPrintf(PETSC_COMM_WORLD, "Setting %d\n", k);CHKERRQ(ierr);
161*4b11644fSPeter Brune     ierr = VecAXPY(x, -1.0, p);CHKERRQ(ierr);
162*4b11644fSPeter Brune   }
163*4b11644fSPeter Brune   if (i == snes->max_its) {
164*4b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
165*4b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
166*4b11644fSPeter Brune   }
167*4b11644fSPeter Brune   PetscFunctionReturn(0);
168*4b11644fSPeter Brune }
169*4b11644fSPeter Brune 
170*4b11644fSPeter Brune 
171*4b11644fSPeter Brune #undef __FUNCT__
172*4b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
173*4b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
174*4b11644fSPeter Brune {
175*4b11644fSPeter Brune   QNContext * qn = (QNContext *)snes->data;
176*4b11644fSPeter Brune   PetscErrorCode ierr;
177*4b11644fSPeter Brune   PetscFunctionBegin;
178*4b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
179*4b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
180*4b11644fSPeter Brune   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
181*4b11644fSPeter Brune   PetscFunctionReturn(0);
182*4b11644fSPeter Brune }
183*4b11644fSPeter Brune 
184*4b11644fSPeter Brune #undef __FUNCT__
185*4b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
186*4b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
187*4b11644fSPeter Brune {
188*4b11644fSPeter Brune   PetscErrorCode ierr;
189*4b11644fSPeter Brune   QNContext * qn;
190*4b11644fSPeter Brune   PetscFunctionBegin;
191*4b11644fSPeter Brune   if (snes->data) {
192*4b11644fSPeter Brune     qn = (QNContext *)snes->data;
193*4b11644fSPeter Brune     if (qn->dX) {
194*4b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
195*4b11644fSPeter Brune     }
196*4b11644fSPeter Brune     if (qn->dF) {
197*4b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
198*4b11644fSPeter Brune     }
199*4b11644fSPeter Brune   }
200*4b11644fSPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
201*4b11644fSPeter Brune   PetscFunctionReturn(0);
202*4b11644fSPeter Brune }
203*4b11644fSPeter Brune 
204*4b11644fSPeter Brune #undef __FUNCT__
205*4b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
206*4b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
207*4b11644fSPeter Brune {
208*4b11644fSPeter Brune   PetscErrorCode ierr;
209*4b11644fSPeter Brune   PetscFunctionBegin;
210*4b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
211*4b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
212*4b11644fSPeter Brune   PetscFunctionReturn(0);
213*4b11644fSPeter Brune }
214*4b11644fSPeter Brune 
215*4b11644fSPeter Brune #undef __FUNCT__
216*4b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
217*4b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
218*4b11644fSPeter Brune {
219*4b11644fSPeter Brune 
220*4b11644fSPeter Brune   PetscErrorCode ierr;
221*4b11644fSPeter Brune   QNContext * qn;
222*4b11644fSPeter Brune 
223*4b11644fSPeter Brune   PetscFunctionBegin;
224*4b11644fSPeter Brune 
225*4b11644fSPeter Brune   qn = (QNContext *)snes->data;
226*4b11644fSPeter Brune 
227*4b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
228*4b11644fSPeter Brune   ierr = PetscOptionsReal("-snes_qn_lambda", "SOR mixing parameter", "SNES", qn->lambda, &qn->lambda, PETSC_NULL);CHKERRQ(ierr);
229*4b11644fSPeter Brune   ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
230*4b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
231*4b11644fSPeter Brune   PetscFunctionReturn(0);
232*4b11644fSPeter Brune }
233*4b11644fSPeter Brune 
234*4b11644fSPeter Brune /* -------------------------------------------------------------------------- */
235*4b11644fSPeter Brune /*MC
236*4b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
237*4b11644fSPeter Brune 
238*4b11644fSPeter Brune 
239*4b11644fSPeter Brune    Level: beginner
240*4b11644fSPeter Brune 
241*4b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
242*4b11644fSPeter Brune M*/
243*4b11644fSPeter Brune EXTERN_C_BEGIN
244*4b11644fSPeter Brune #undef __FUNCT__
245*4b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
246*4b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
247*4b11644fSPeter Brune {
248*4b11644fSPeter Brune 
249*4b11644fSPeter Brune   PetscErrorCode ierr;
250*4b11644fSPeter Brune   QNContext * qn;
251*4b11644fSPeter Brune 
252*4b11644fSPeter Brune   PetscFunctionBegin;
253*4b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
254*4b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
255*4b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
256*4b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
257*4b11644fSPeter Brune   snes->ops->view            = 0;
258*4b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
259*4b11644fSPeter Brune 
260*4b11644fSPeter Brune   ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr);
261*4b11644fSPeter Brune   snes->data = (void *) qn;
262*4b11644fSPeter Brune   qn->m = 10;
263*4b11644fSPeter Brune   qn->lambda = 1.;
264*4b11644fSPeter Brune   qn->dX = PETSC_NULL;
265*4b11644fSPeter Brune   qn->dF = PETSC_NULL;
266*4b11644fSPeter Brune   PetscFunctionReturn(0);
267*4b11644fSPeter Brune }
268*4b11644fSPeter Brune EXTERN_C_END
269