Actual source code: sorqn.c
2: #include <private/snesimpl.h>
4: typedef struct {
5: PetscBool jacobian_start; /* start with Bi = Jii */
6: PetscInt n_restart; /* after n iterations, set Bi_n = Jii(x)_n */
7: PetscReal alpha; /* SOR mixing parameter */
8: } SORQNContext;
12: static PetscErrorCode SNESSolve_SORQN(SNES snes)
13: {
14: PetscErrorCode ierr;
15: SORQNContext * sorqn;
16: PetscReal f_norm, f_norm_old;
17: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
18: Vec X, F, dX, Y, B;
19: PetscInt i, j;
20: PetscInt rs, re;
21: PetscScalar dX_i, Y_i;
24: snes->reason = SNES_CONVERGED_ITERATING;
25: sorqn = (SORQNContext *)snes->data;
27: X = snes->vec_sol;
28: F = snes->vec_func;
29: Y = snes->work[1];
30: dX = snes->vec_sol_update;
32: /* the diagonal of the approximate jacobian (broyden update) */
33: B = snes->work[0];
35: VecGetOwnershipRange(X, &rs, &re);
37: /* initial residual */
38: SNESComputeFunction(snes,X,F);
39: if (snes->domainerror) {
40: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
41: return(0);
42: }
44: /* set the initial guess for the broyden update */
45: if (sorqn->jacobian_start) {
46: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
47: MatGetDiagonal(snes->jacobian, B);
48: } else {
49: VecSet(B, 1.0);
50: }
51:
52: /* take the initial residual F and residual norm*/
53: SNESComputeFunction(snes, X, F);
54: VecNorm(F, NORM_2, &f_norm);
55: PetscObjectTakeAccess(snes);
56: snes->norm = f_norm;
57: snes->iter = 0;
58: PetscObjectGrantAccess(snes);
59: SNESLogConvHistory(snes, f_norm, 0);
60: SNESMonitor(snes, 0, f_norm);
61: (*snes->ops->converged)(snes,0,0.0,0.0,f_norm,&snes->reason,snes->cnvP);
62: if (snes->reason) return(0);
63: for (i = 1; i < snes->max_its; i++) {
64: f_norm_old = f_norm;
65: /* Compute the change in X */
66: VecCopy(F, dX);
67: VecScale(dX, sorqn->alpha);
68: VecPointwiseDivide(dX, dX, B);
69: VecAXPY(X, -1.0, dX);
70:
71: /* Compute the update for B */
72: VecCopy(F, Y);
73: SNESComputeFunction(snes, X, F); /* r^{i+1} = F(x^i) */
74: VecNorm(F, NORM_2, &f_norm);
75: /*check the convergence */
76: PetscObjectTakeAccess(snes);
77: snes->iter = i;
78: snes->norm = f_norm;
79: PetscObjectGrantAccess(snes);
80: SNESLogConvHistory(snes, f_norm, i);
81: SNESMonitor(snes, i, f_norm);
82: (*snes->ops->converged)(snes,i,0.0,0.0,f_norm,&snes->reason,snes->cnvP);
83: if (snes->reason) {
84: return(0);
85: }
86: VecAXPBY(Y, -1.0, 1.0, F); /* y^i = r^{i} - r^{i-1} */
87:
88: /*CHOICE: either restart, or continue doing diagonal rank-one updates */
89:
90: if (i % sorqn->n_restart == 0 || f_norm > 2.0*f_norm_old) {
91: if (sorqn->jacobian_start) {
92: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
93: MatGetDiagonal(snes->jacobian, B);
94: } else {
95: VecSet(B, 1.0);
96: }
97: } else {
98: for (j = rs; j < re; j++) {
99: VecGetValues(dX, 1, &j, &dX_i);
100: VecGetValues(Y, 1, &j, &Y_i);
101: #ifdef PETSC_USE_COMPLEX
102: if (PetscAbs(PetscRealPart(dX_i)) > 1e-18) {
103: Y_i = Y_i / dX_i;
104: if (PetscAbs(PetscRealPart(Y_i)) > 1e-6) {
105: VecSetValues(B, 1, &j, &Y_i, INSERT_VALUES);
106: }
107: }
108: #else
109: if (PetscAbs(dX_i) > 1e-18) {
110: Y_i = Y_i / dX_i;
111: if (PetscAbs(Y_i) > 1e-6) {
112: VecSetValues(B, 1, &j, &Y_i, INSERT_VALUES);
113: }
114: }
115: #endif
116: }
117: VecAssemblyBegin(B);
118: VecAssemblyEnd(B);
119: }
120: }
121: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
122: return(0);
123: }
127: static PetscErrorCode SNESSetUp_SORQN(SNES snes)
128: {
131: SNESDefaultGetWork(snes,2);
132: return(0);
133: }
137: static PetscErrorCode SNESReset_SORQN(SNES snes)
138: {
141: if (snes->work) {VecDestroyVecs(snes->nwork,&snes->work);}
142: return(0);
143: }
147: static PetscErrorCode SNESDestroy_SORQN(SNES snes)
148: {
151: SNESReset_SORQN(snes);
152: PetscFree(snes->data);
153: return(0);
154: }
158: static PetscErrorCode SNESSetFromOptions_SORQN(SNES snes)
159: {
162: SORQNContext * sorqn;
166: sorqn = (SORQNContext *)snes->data;
168: PetscOptionsHead("SNES SOR-QN options");
169: PetscOptionsBool("-snes_sorqn_jacobian_start", "Start Quasi-Newton with actual Jacobian", "SNES", sorqn->jacobian_start, &sorqn->jacobian_start, PETSC_NULL);
170: PetscOptionsReal("-snes_sorqn_alpha", "SOR mixing parameter", "SNES", sorqn->alpha, &sorqn->alpha, PETSC_NULL);
171: PetscOptionsInt("-snes_sorqn_restart", "Iterations before Newton restart", "SNES", sorqn->n_restart, &sorqn->n_restart, PETSC_NULL);
172: PetscOptionsTail();
173: return(0);
174: }
176: /* -------------------------------------------------------------------------- */
177: /*MC
178: SNESSORQN - Nonlinear solver based upon the SOR-Quasi-Newton method.
179: Reference:
181: Martinez, J.M., SOR-Secant methods,
182: SIAM Journal on Numerical Analysis, Vol. 31, No. 1 (Feb. 1994), SIAM
184: This implementation is still very experimental, and needs to be modified to use
185: the inner quasi-Newton iteration on blocks of unknowns.
187: Level: beginner
189: .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
190: M*/
194: PetscErrorCode SNESCreate_SORQN(SNES snes)
195: {
196:
198: SORQNContext * sorqn;
201: snes->ops->setup = SNESSetUp_SORQN;
202: snes->ops->solve = SNESSolve_SORQN;
203: snes->ops->destroy = SNESDestroy_SORQN;
204: snes->ops->setfromoptions = SNESSetFromOptions_SORQN;
205: snes->ops->view = 0;
206: snes->ops->reset = SNESReset_SORQN;
208: PetscNewLog(snes, SORQNContext, &sorqn);
209: snes->data = (void *) sorqn;
210: sorqn->jacobian_start = PETSC_FALSE;
211: sorqn->alpha = 0.25;
212: sorqn->n_restart = 10;
213: return(0);
214: }