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