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