xref: /petsc/src/snes/impls/gs/snesgs.c (revision f692024ea6ceda132bc9ff30ca7a31e85cfbbcb2)
1 #include <private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2 
3 typedef struct {
4   PetscInt sweeps;
5 } SNES_GS;
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "SNESReset_GS"
9 PetscErrorCode SNESReset_GS(SNES snes)
10 {
11   PetscFunctionBegin;
12   PetscFunctionReturn(0);
13 }
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "SNESDestroy_GS"
17 PetscErrorCode SNESDestroy_GS(SNES snes)
18 {
19   PetscErrorCode ierr;
20 
21   PetscFunctionBegin;
22   ierr = SNESReset_GS(snes);CHKERRQ(ierr);
23   ierr = PetscFree(snes->data);
24   PetscFunctionReturn(0);
25 }
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "SNESSetUp_GS"
29 PetscErrorCode SNESSetUp_GS(SNES snes)
30 {
31   PetscFunctionBegin;
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "SNESSetFromOptions_GS"
37 PetscErrorCode SNESSetFromOptions_GS(SNES snes)
38 {
39   PetscErrorCode ierr;
40 
41   PetscFunctionBegin;
42   ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr);
43   PetscFunctionReturn(0);
44 }
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "SNESView_GS"
48 PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer)
49 {
50   PetscFunctionBegin;
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "SNESSolve_GS"
56 PetscErrorCode SNESSolve_GS(SNES snes)
57 {
58   Vec            F;
59   Vec            X;
60   Vec            B;
61   PetscInt       i;
62   PetscReal      fnorm;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66 
67   X = snes->vec_sol;
68   F = snes->vec_func;
69   B = snes->vec_rhs;
70 
71   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
72   snes->iter = 0;
73   snes->norm = 0.;
74   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
75 
76   snes->reason = SNES_CONVERGED_ITS;
77   /* compute the initial function and preconditioned update delX */
78   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
79   if (snes->domainerror) {
80     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
81     PetscFunctionReturn(0);
82   }
83   /* convergence test */
84   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
85   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
86   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
87   snes->norm = fnorm;
88   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
89   SNESLogConvHistory(snes,fnorm,0);
90   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
91 
92   /* set parameter for default relative tolerance convergence test */
93   snes->ttol = fnorm*snes->rtol;
94   /* test convergence */
95   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
96   if (snes->reason) PetscFunctionReturn(0);
97 
98   /* Call general purpose update function */
99   if (snes->ops->update) {
100     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
101  }
102   for(i = 0; i < snes->max_its; i++) {
103     ierr = SNESComputeGS(snes, B, X);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     /* Monitor convergence */
112     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
113     snes->iter = i+1;
114     snes->norm = fnorm;
115     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
116     SNESLogConvHistory(snes,snes->norm,0);
117     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
118 
119     /* Test for convergence */
120     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
121     if (snes->reason) PetscFunctionReturn(0);
122 
123     /* Call general purpose update function */
124     if (snes->ops->update) {
125       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
126     }
127   }
128   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
129   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
130   PetscFunctionReturn(0);
131 }
132 
133 /*MC
134   SNESGS - Just calls the user-provided solution routine provided with SNESSetGS()
135 
136    Level: advanced
137 
138   Notes:
139   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetPC(), it will have
140   its parent's Gauss-Seidel routine associated with it.
141 
142 .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types)
143 M*/
144 
145 EXTERN_C_BEGIN
146 #undef __FUNCT__
147 #define __FUNCT__ "SNESCreate_GS"
148 PetscErrorCode SNESCreate_GS(SNES snes)
149 {
150   SNES_GS        *gs;
151   PetscErrorCode ierr;
152 
153   PetscFunctionBegin;
154   snes->ops->destroy        = SNESDestroy_GS;
155   snes->ops->setup          = SNESSetUp_GS;
156   snes->ops->setfromoptions = SNESSetFromOptions_GS;
157   snes->ops->view           = SNESView_GS;
158   snes->ops->solve          = SNESSolve_GS;
159   snes->ops->reset          = SNESReset_GS;
160 
161   snes->usesksp             = PETSC_FALSE;
162   snes->usespc              = PETSC_FALSE;
163 
164   ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr);
165   snes->data = (void*) gs;
166   PetscFunctionReturn(0);
167 }
168 EXTERN_C_END
169