xref: /petsc/src/snes/impls/gs/snesgs.c (revision 6cd56b8b215a5e64de6d370e6130585cb46bf7a3)
1 #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2 
3 typedef struct {
4   PetscInt  sweeps;
5 } SNES_GS;
6 
7 
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "SNESGSSetSweeps"
11 /*@
12    SNESGSSetSweeps - Sets the number of sweeps of GS to use.
13 
14    Input Parameters:
15 +  snes   - the SNES context
16 -  sweeps  - the number of sweeps of GS to perform.
17 
18    Level: intermediate
19 
20 .keywords: SNES, nonlinear, set, Gauss-Siedel
21 
22 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps()
23 @*/
24 
25 PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps)
26 {
27   SNES_GS *gs = (SNES_GS*)snes->data;
28 
29   PetscFunctionBegin;
30   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
31   gs->sweeps = sweeps;
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "SNESGSGetSweeps"
37 /*@
38    SNESGSGetSweeps - Gets the number of sweeps GS will use.
39 
40    Input Parameters:
41 .  snes   - the SNES context
42 
43    Output Parameters:
44 .  sweeps  - the number of sweeps of GS to perform.
45 
46    Level: intermediate
47 
48 .keywords: SNES, nonlinear, set, Gauss-Siedel
49 
50 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps()
51 @*/
52 PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps)
53 {
54   SNES_GS *gs = (SNES_GS*)snes->data;
55 
56   PetscFunctionBegin;
57   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
58   *sweeps = gs->sweeps;
59   PetscFunctionReturn(0);
60 }
61 
62 
63 #undef __FUNCT__
64 #define __FUNCT__ "SNESDefaultApplyGS"
65 PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx)
66 {
67   PetscFunctionBegin;
68   /* see if there's a coloring on the DM */
69 
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "SNESReset_GS"
75 PetscErrorCode SNESReset_GS(SNES snes)
76 {
77   PetscFunctionBegin;
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "SNESDestroy_GS"
83 PetscErrorCode SNESDestroy_GS(SNES snes)
84 {
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88   ierr = SNESReset_GS(snes);CHKERRQ(ierr);
89   ierr = PetscFree(snes->data);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "SNESSetUp_GS"
95 PetscErrorCode SNESSetUp_GS(SNES snes)
96 {
97   PetscFunctionBegin;
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "SNESSetFromOptions_GS"
103 PetscErrorCode SNESSetFromOptions_GS(SNES snes)
104 {
105   SNES_GS          *gs = (SNES_GS*)snes->data;
106   PetscErrorCode   ierr;
107 
108   PetscFunctionBegin;
109   ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr);
110   /* GS Options */
111   ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&gs->sweeps,PETSC_NULL);CHKERRQ(ierr);
112   ierr = PetscOptionsTail();CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "SNESView_GS"
118 PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer)
119 {
120   PetscFunctionBegin;
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "SNESSolve_GS"
126 PetscErrorCode SNESSolve_GS(SNES snes)
127 {
128   Vec            F;
129   Vec            X;
130   Vec            B;
131   PetscInt       i;
132   PetscReal      fnorm;
133   PetscErrorCode ierr;
134   SNESNormType   normtype;
135 
136   PetscFunctionBegin;
137 
138   X = snes->vec_sol;
139   F = snes->vec_func;
140   B = snes->vec_rhs;
141 
142   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
143   snes->iter = 0;
144   snes->norm = 0.;
145   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
146   snes->reason = SNES_CONVERGED_ITERATING;
147 
148   ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
149   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
150     /* compute the initial function and preconditioned update delX */
151     if (!snes->vec_func_init_set) {
152       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
153       if (snes->domainerror) {
154         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
155         PetscFunctionReturn(0);
156       }
157     } else {
158       snes->vec_func_init_set = PETSC_FALSE;
159     }
160 
161     /* convergence test */
162     if (!snes->norm_init_set) {
163       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
164       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
165     } else {
166       fnorm = snes->norm_init;
167       snes->norm_init_set = PETSC_FALSE;
168     }
169     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
170     snes->iter = 0;
171     snes->norm = fnorm;
172     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
173     SNESLogConvHistory(snes,snes->norm,0);
174     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
175 
176     /* set parameter for default relative tolerance convergence test */
177     snes->ttol = fnorm*snes->rtol;
178 
179     /* test convergence */
180     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
181     if (snes->reason) PetscFunctionReturn(0);
182   } else {
183     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
184     SNESLogConvHistory(snes,snes->norm,0);
185     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
186   }
187 
188 
189   /* Call general purpose update function */
190   if (snes->ops->update) {
191     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
192   }
193 
194   for (i = 0; i < snes->max_its; i++) {
195     ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
196     /* only compute norms if requested or about to exit due to maximum iterations */
197     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
198       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
199       if (snes->domainerror) {
200         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
201         PetscFunctionReturn(0);
202       }
203       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
204       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
205     }
206     /* Monitor convergence */
207     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
208     snes->iter = i+1;
209     snes->norm = fnorm;
210     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
211     SNESLogConvHistory(snes,snes->norm,0);
212     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
213     /* Test for convergence */
214     if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
215     if (snes->reason) PetscFunctionReturn(0);
216     /* Call general purpose update function */
217     if (snes->ops->update) {
218       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
219     }
220   }
221   if (normtype == SNES_NORM_FUNCTION) {
222     if (i == snes->max_its) {
223       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
224       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
225     }
226   } else {
227     if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
228   }
229   PetscFunctionReturn(0);
230 }
231 
232 /*MC
233   SNESGS - Just calls the user-provided solution routine provided with SNESSetGS()
234 
235    Level: advanced
236 
237   Notes:
238   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetPC(), it will have
239   its parent's Gauss-Seidel routine associated with it.
240 
241 .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types)
242 M*/
243 
244 EXTERN_C_BEGIN
245 #undef __FUNCT__
246 #define __FUNCT__ "SNESCreate_GS"
247 PetscErrorCode SNESCreate_GS(SNES snes)
248 {
249   SNES_GS        *gs;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   snes->ops->destroy        = SNESDestroy_GS;
254   snes->ops->setup          = SNESSetUp_GS;
255   snes->ops->setfromoptions = SNESSetFromOptions_GS;
256   snes->ops->view           = SNESView_GS;
257   snes->ops->solve          = SNESSolve_GS;
258   snes->ops->reset          = SNESReset_GS;
259 
260   snes->usesksp             = PETSC_FALSE;
261   snes->usespc              = PETSC_FALSE;
262 
263   if (!snes->tolerancesset) {
264     snes->max_its             = 10000;
265     snes->max_funcs           = 10000;
266   }
267 
268   ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr);
269 
270   gs->sweeps = 1;
271 
272   snes->data = (void*) gs;
273 
274   PetscFunctionReturn(0);
275 }
276 EXTERN_C_END
277