xref: /petsc/src/snes/impls/composite/snescomposite.c (revision 90a8ba9be28a04c4922db9ec087f82b84f0205dc)
1eed5f15bSPeter Brune 
2eed5f15bSPeter Brune /*
3eed5f15bSPeter Brune       Defines a SNES that can consist of a collection of SNESes
4eed5f15bSPeter Brune */
5eed5f15bSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
6*90a8ba9bSPeter Brune #include <petscblaslapack.h>
7eed5f15bSPeter Brune 
8*90a8ba9bSPeter Brune const char *const        SNESCompositeTypes[]   = {"ADDITIVE","MULTIPLICATIVE","ADDITIVEOPTIMAL","SNESCompositeType","SNES_COMPOSITE",0};
9eed5f15bSPeter Brune 
10eed5f15bSPeter Brune typedef struct _SNES_CompositeLink *SNES_CompositeLink;
11eed5f15bSPeter Brune struct _SNES_CompositeLink {
12eed5f15bSPeter Brune   SNES               snes;
138f626970SPeter Brune   PetscReal          dmp;
14*90a8ba9bSPeter Brune   Vec                X;
15eed5f15bSPeter Brune   SNES_CompositeLink next;
16eed5f15bSPeter Brune   SNES_CompositeLink previous;
17eed5f15bSPeter Brune };
18eed5f15bSPeter Brune 
19eed5f15bSPeter Brune typedef struct {
20eed5f15bSPeter Brune   SNES_CompositeLink head;
21*90a8ba9bSPeter Brune   PetscInt           nsnes;
22eed5f15bSPeter Brune   SNESCompositeType  type;
23eed5f15bSPeter Brune   Vec                Xorig;
24*90a8ba9bSPeter Brune 
25*90a8ba9bSPeter Brune   /* context for ADDITIVEOPTIMAL */
26*90a8ba9bSPeter Brune   Vec                *Xes,*Fes;      /* solution and residual vectors for the subsolvers */
27*90a8ba9bSPeter Brune   PetscScalar        *h;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
28*90a8ba9bSPeter Brune   PetscBLASInt       m;              /* matrix dimension -- nsnes */
29*90a8ba9bSPeter Brune   PetscBLASInt       n;              /* matrix dimension -- nsnes + 1 */
30*90a8ba9bSPeter Brune   PetscBLASInt       nrhs;           /* the number of right hand sides */
31*90a8ba9bSPeter Brune   PetscBLASInt       lda;            /* the padded matrix dimension */
32*90a8ba9bSPeter Brune   PetscBLASInt       ldb;            /* the padded vector dimension */
33*90a8ba9bSPeter Brune   PetscReal          *s;             /* the singular values */
34*90a8ba9bSPeter Brune   PetscReal          *beta;          /* the RHS and combination */
35*90a8ba9bSPeter Brune   PetscReal          rcond;          /* the exit condition */
36*90a8ba9bSPeter Brune   PetscBLASInt       rank;           /* the effective rank */
37*90a8ba9bSPeter Brune   PetscScalar        *work;          /* the work vector */
38*90a8ba9bSPeter Brune   PetscReal          *rwork;         /* the real work vector used for complex */
39*90a8ba9bSPeter Brune   PetscBLASInt       lwork;          /* the size of the work vector */
40*90a8ba9bSPeter Brune   PetscBLASInt       info;           /* the output condition */
41*90a8ba9bSPeter Brune 
42eed5f15bSPeter Brune } SNES_Composite;
43eed5f15bSPeter Brune 
44eed5f15bSPeter Brune #undef __FUNCT__
45eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative"
46eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
47eed5f15bSPeter Brune {
48eed5f15bSPeter Brune   PetscErrorCode     ierr;
49eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
50eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
51eed5f15bSPeter Brune   Vec                FSub;
52eed5f15bSPeter Brune   PetscReal          fsubnorm;
53eed5f15bSPeter Brune 
54eed5f15bSPeter Brune   PetscFunctionBegin;
55eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
56eed5f15bSPeter Brune   if (snes->normtype == SNES_NORM_FUNCTION) {
57eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
58eed5f15bSPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
59eed5f15bSPeter Brune   }
60eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
61eed5f15bSPeter Brune 
62eed5f15bSPeter Brune   while (next->next) {
63eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
64eed5f15bSPeter Brune     if (next->snes->pcside == PC_RIGHT && next->snes->normtype != SNES_NORM_NONE) {
65eed5f15bSPeter Brune       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
66eed5f15bSPeter Brune       ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr);
67eed5f15bSPeter Brune       next = next->next;
68eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
69eed5f15bSPeter Brune       ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr);
70eed5f15bSPeter Brune     } else {
71eed5f15bSPeter Brune       next = next->next;
72eed5f15bSPeter Brune     }
73eed5f15bSPeter Brune     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
74eed5f15bSPeter Brune   }
75eed5f15bSPeter Brune   if (next->snes->pcside == PC_RIGHT) {
76eed5f15bSPeter Brune     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
77eed5f15bSPeter Brune     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
78eed5f15bSPeter Brune     if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);}
79eed5f15bSPeter Brune   } else if (snes->normtype == SNES_NORM_FUNCTION) {
80eed5f15bSPeter Brune     SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
81eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
82eed5f15bSPeter Brune   }
83eed5f15bSPeter Brune   PetscFunctionReturn(0);
84eed5f15bSPeter Brune }
85eed5f15bSPeter Brune 
86eed5f15bSPeter Brune #undef __FUNCT__
87eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive"
88eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
89eed5f15bSPeter Brune {
90eed5f15bSPeter Brune   PetscErrorCode     ierr;
91eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
92eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
93eed5f15bSPeter Brune   Vec                Y,Xorig;
94eed5f15bSPeter Brune 
95eed5f15bSPeter Brune   PetscFunctionBegin;
96eed5f15bSPeter Brune   Y = snes->vec_sol_update;
97eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
98eed5f15bSPeter Brune   Xorig = jac->Xorig;
99eed5f15bSPeter Brune   ierr = VecCopy(X,Xorig);
100eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
101eed5f15bSPeter Brune   if (snes->normtype == SNES_NORM_FUNCTION) {
102eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
103eed5f15bSPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
104eed5f15bSPeter Brune     while (next->next) {
105eed5f15bSPeter Brune       next = next->next;
106eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
107eed5f15bSPeter Brune       if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
108eed5f15bSPeter Brune     }
109eed5f15bSPeter Brune   }
110eed5f15bSPeter Brune   next = jac->head;
1118f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
112eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
113eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1148f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1158f626970SPeter Brune   while (next->next) {
1168f626970SPeter Brune     next = next->next;
1178f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1188f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
1198f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1208f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
121eed5f15bSPeter Brune   }
122eed5f15bSPeter Brune   if (snes->normtype == SNES_NORM_FUNCTION) {
123eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
124eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
125eed5f15bSPeter Brune   }
126eed5f15bSPeter Brune   PetscFunctionReturn(0);
127eed5f15bSPeter Brune }
128*90a8ba9bSPeter Brune 
129*90a8ba9bSPeter Brune #undef __FUNCT__
130*90a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal"
131*90a8ba9bSPeter Brune /* approximately solve the overdetermined system:
132*90a8ba9bSPeter Brune 
133*90a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
134*90a8ba9bSPeter Brune  \alpha_i                      = 1
135*90a8ba9bSPeter Brune 
136*90a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
137*90a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
138*90a8ba9bSPeter Brune 
139*90a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
140*90a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
141*90a8ba9bSPeter Brune  */
142*90a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
143*90a8ba9bSPeter Brune {
144*90a8ba9bSPeter Brune   PetscErrorCode     ierr;
145*90a8ba9bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
146*90a8ba9bSPeter Brune   SNES_CompositeLink next = jac->head;
147*90a8ba9bSPeter Brune   Vec                *Xes = jac->Xes,*Fes = jac->Fes;
148*90a8ba9bSPeter Brune   PetscInt           i,j;
149*90a8ba9bSPeter Brune   PetscScalar        maxdiag = 0.0;
150*90a8ba9bSPeter Brune 
151*90a8ba9bSPeter Brune   PetscFunctionBegin;
152*90a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
153*90a8ba9bSPeter Brune   next = jac->head;
154*90a8ba9bSPeter Brune   i = 0;
155*90a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
156*90a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
157*90a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,Xes[i],Fes[i]);CHKERRQ(ierr);
158*90a8ba9bSPeter Brune   while (next->next) {
159*90a8ba9bSPeter Brune     i++;
160*90a8ba9bSPeter Brune     next = next->next;
161*90a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
162*90a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
163*90a8ba9bSPeter Brune     ierr = SNESComputeFunction(snes,Xes[i],Fes[i]);CHKERRQ(ierr);
164*90a8ba9bSPeter Brune   }
165*90a8ba9bSPeter Brune 
166*90a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
167*90a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
168*90a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
169*90a8ba9bSPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->m]);CHKERRQ(ierr);
170*90a8ba9bSPeter Brune     }
171*90a8ba9bSPeter Brune   }
172*90a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
173*90a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
174*90a8ba9bSPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->m]);CHKERRQ(ierr);
175*90a8ba9bSPeter Brune       jac->h[i + j*jac->m] *= 2.;
176*90a8ba9bSPeter Brune     }
177*90a8ba9bSPeter Brune   }
178*90a8ba9bSPeter Brune 
179*90a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
180*90a8ba9bSPeter Brune     if (maxdiag < jac->h[i+i*jac->m] && jac->h[i+i*jac->m] > 0.) maxdiag = jac->h[i+i*jac->m];
181*90a8ba9bSPeter Brune   }
182*90a8ba9bSPeter Brune 
183*90a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
184*90a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
185*90a8ba9bSPeter Brune       jac->h[i + j*jac->m] = jac->h[j + i*jac->m];
186*90a8ba9bSPeter Brune     }
187*90a8ba9bSPeter Brune   }
188*90a8ba9bSPeter Brune 
189*90a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
190*90a8ba9bSPeter Brune     jac->beta[i] = 0.;
191*90a8ba9bSPeter Brune     jac->h[jac->m-1 + jac->m*i] = 2.*maxdiag;
192*90a8ba9bSPeter Brune   }
193*90a8ba9bSPeter Brune   jac->beta[jac->m-1] = 2.*maxdiag;
194*90a8ba9bSPeter Brune 
195*90a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
196*90a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
197*90a8ba9bSPeter Brune #else
198*90a8ba9bSPeter Brune   jac->info  = 0;
199*90a8ba9bSPeter Brune   jac->rcond = -1.;
200*90a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
201*90a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
202*90a8ba9bSPeter Brune   PetscStackCall("LAPACKgelss",LAPACKgelss_(&jac->m,&jac->n,&jac->nrhs,jac->h,&jac->lda,jac->beta,&jac->lda,jac->s,&jac->rcond,&jac->rank,jac->work,&jac->lwork,jac->rwork,&jac->info));
203*90a8ba9bSPeter Brune #else
204*90a8ba9bSPeter Brune   PetscStackCall("LAPACKgelss",LAPACKgelss_(&jac->m,&jac->n,&jac->nrhs,jac->h,&jac->lda,jac->beta,&jac->lda,jac->s,&jac->rcond,&jac->rank,jac->work,&jac->lwork,&jac->info));
205*90a8ba9bSPeter Brune #endif
206*90a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
207*90a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
208*90a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
209*90a8ba9bSPeter Brune #endif
210*90a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
211*90a8ba9bSPeter Brune     if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
212*90a8ba9bSPeter Brune   }
213*90a8ba9bSPeter Brune   ierr = VecSet(X,0.);CHKERRQ(ierr);
214*90a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
215*90a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
216*90a8ba9bSPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
217*90a8ba9bSPeter Brune 
218*90a8ba9bSPeter Brune   PetscFunctionReturn(0);
219*90a8ba9bSPeter Brune }
220*90a8ba9bSPeter Brune 
221eed5f15bSPeter Brune #undef __FUNCT__
222eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite"
223eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
224eed5f15bSPeter Brune {
225eed5f15bSPeter Brune   PetscErrorCode     ierr;
226eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
227eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
228*90a8ba9bSPeter Brune   PetscInt           n=0,i;
229*90a8ba9bSPeter Brune   Vec                F;
230eed5f15bSPeter Brune 
231eed5f15bSPeter Brune   PetscFunctionBegin;
232eed5f15bSPeter Brune   while (next) {
233*90a8ba9bSPeter Brune     n++;
234eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
235eed5f15bSPeter Brune     next = next->next;
236eed5f15bSPeter Brune   }
237*90a8ba9bSPeter Brune   jac->nsnes = n;
238*90a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
239*90a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
240*90a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
241*90a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr);
242*90a8ba9bSPeter Brune     next = jac->head;
243*90a8ba9bSPeter Brune     i = 0;
244*90a8ba9bSPeter Brune     while (next) {
245*90a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
246*90a8ba9bSPeter Brune       jac->Fes[i] = F;
247*90a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
248*90a8ba9bSPeter Brune       next = next->next;
249*90a8ba9bSPeter Brune       i++;
250*90a8ba9bSPeter Brune     }
251*90a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
252*90a8ba9bSPeter Brune     jac->nrhs  = 1;
253*90a8ba9bSPeter Brune     jac->lda   = jac->nsnes+1;
254*90a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
255*90a8ba9bSPeter Brune     jac->m     = jac->nsnes+1;
256*90a8ba9bSPeter Brune     jac->n     = jac->nsnes;
257*90a8ba9bSPeter Brune 
258*90a8ba9bSPeter Brune     ierr = PetscMalloc(jac->m*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr);
259*90a8ba9bSPeter Brune     ierr = PetscMalloc(jac->m*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr);
260*90a8ba9bSPeter Brune     ierr = PetscMalloc(jac->m*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr);
261*90a8ba9bSPeter Brune     jac->lwork = 12*jac->m;
262*90a8ba9bSPeter Brune #if PETSC_USE_COMPLEX
263*90a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr);
264*90a8ba9bSPeter Brune #endif
265*90a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr);
266*90a8ba9bSPeter Brune   }
267*90a8ba9bSPeter Brune 
268eed5f15bSPeter Brune   PetscFunctionReturn(0);
269eed5f15bSPeter Brune }
270eed5f15bSPeter Brune 
271eed5f15bSPeter Brune #undef __FUNCT__
272eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite"
273eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
274eed5f15bSPeter Brune {
275eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
276eed5f15bSPeter Brune   PetscErrorCode   ierr;
277eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
278eed5f15bSPeter Brune 
279eed5f15bSPeter Brune   PetscFunctionBegin;
280eed5f15bSPeter Brune   while (next) {
281eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
282eed5f15bSPeter Brune     next = next->next;
283eed5f15bSPeter Brune   }
284eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
285*90a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
286*90a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
287*90a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
288*90a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
289*90a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
290*90a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
291*90a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
292*90a8ba9bSPeter Brune 
293eed5f15bSPeter Brune   PetscFunctionReturn(0);
294eed5f15bSPeter Brune }
295eed5f15bSPeter Brune 
296eed5f15bSPeter Brune #undef __FUNCT__
297eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite"
298eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
299eed5f15bSPeter Brune {
300eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
301eed5f15bSPeter Brune   PetscErrorCode   ierr;
302eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
303eed5f15bSPeter Brune 
304eed5f15bSPeter Brune   PetscFunctionBegin;
305eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
306eed5f15bSPeter Brune   while (next) {
307eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
308eed5f15bSPeter Brune     next_tmp = next;
309eed5f15bSPeter Brune     next     = next->next;
310eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
311eed5f15bSPeter Brune   }
312eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
313eed5f15bSPeter Brune   PetscFunctionReturn(0);
314eed5f15bSPeter Brune }
315eed5f15bSPeter Brune 
316eed5f15bSPeter Brune #undef __FUNCT__
317eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite"
318eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
319eed5f15bSPeter Brune {
320eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
321eed5f15bSPeter Brune   PetscErrorCode     ierr;
322eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
323eed5f15bSPeter Brune   SNES_CompositeLink next;
324eed5f15bSPeter Brune   char               *sneses[8];
3258f626970SPeter Brune   PetscReal          dmps[8];
326eed5f15bSPeter Brune   PetscBool          flg;
327eed5f15bSPeter Brune 
328eed5f15bSPeter Brune   PetscFunctionBegin;
329eed5f15bSPeter Brune   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
330eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
331eed5f15bSPeter Brune   if (flg) {
332eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
333eed5f15bSPeter Brune   }
334eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
335eed5f15bSPeter Brune   if (flg) {
336eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
337eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
338eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
339eed5f15bSPeter Brune     }
340eed5f15bSPeter Brune   }
3418f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
3428f626970SPeter Brune   if (flg) {
3438f626970SPeter Brune     for (i=0; i<nmax; i++) {
3448f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
3458f626970SPeter Brune     }
3468f626970SPeter Brune   }
347eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
348eed5f15bSPeter Brune 
349eed5f15bSPeter Brune   next = jac->head;
350eed5f15bSPeter Brune   while (next) {
351eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
352eed5f15bSPeter Brune     next = next->next;
353eed5f15bSPeter Brune   }
354eed5f15bSPeter Brune   PetscFunctionReturn(0);
355eed5f15bSPeter Brune }
356eed5f15bSPeter Brune 
357eed5f15bSPeter Brune #undef __FUNCT__
358eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite"
359eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
360eed5f15bSPeter Brune {
361eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
362eed5f15bSPeter Brune   PetscErrorCode   ierr;
363eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
364eed5f15bSPeter Brune   PetscBool        iascii;
365eed5f15bSPeter Brune 
366eed5f15bSPeter Brune   PetscFunctionBegin;
367eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
368eed5f15bSPeter Brune   if (iascii) {
369eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
370eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
371eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
372eed5f15bSPeter Brune   }
373eed5f15bSPeter Brune   if (iascii) {
374eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
375eed5f15bSPeter Brune   }
376eed5f15bSPeter Brune   while (next) {
377eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
378eed5f15bSPeter Brune     next = next->next;
379eed5f15bSPeter Brune   }
380eed5f15bSPeter Brune   if (iascii) {
381eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
382eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
383eed5f15bSPeter Brune   }
384eed5f15bSPeter Brune   PetscFunctionReturn(0);
385eed5f15bSPeter Brune }
386eed5f15bSPeter Brune 
387eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
388eed5f15bSPeter Brune 
389eed5f15bSPeter Brune #undef __FUNCT__
390eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite"
391eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
392eed5f15bSPeter Brune {
393eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
394eed5f15bSPeter Brune 
395eed5f15bSPeter Brune   PetscFunctionBegin;
396eed5f15bSPeter Brune   jac->type = type;
397eed5f15bSPeter Brune   PetscFunctionReturn(0);
398eed5f15bSPeter Brune }
399eed5f15bSPeter Brune 
400eed5f15bSPeter Brune #undef __FUNCT__
401eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite"
402eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
403eed5f15bSPeter Brune {
404eed5f15bSPeter Brune   SNES_Composite     *jac;
405eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
406eed5f15bSPeter Brune   PetscErrorCode   ierr;
407eed5f15bSPeter Brune   PetscInt         cnt = 0;
408eed5f15bSPeter Brune   const char       *prefix;
409eed5f15bSPeter Brune   char             newprefix[8];
410eed5f15bSPeter Brune   DM               dm;
411eed5f15bSPeter Brune 
412eed5f15bSPeter Brune   PetscFunctionBegin;
413eed5f15bSPeter Brune   ierr        = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr);
414eed5f15bSPeter Brune   ilink->next = 0;
415eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
416eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
417eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
418eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
419eed5f15bSPeter Brune 
420eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
421eed5f15bSPeter Brune   next = jac->head;
422eed5f15bSPeter Brune   if (!next) {
423eed5f15bSPeter Brune     jac->head       = ilink;
424eed5f15bSPeter Brune     ilink->previous = NULL;
425eed5f15bSPeter Brune   } else {
426eed5f15bSPeter Brune     cnt++;
427eed5f15bSPeter Brune     while (next->next) {
428eed5f15bSPeter Brune       next = next->next;
429eed5f15bSPeter Brune       cnt++;
430eed5f15bSPeter Brune     }
431eed5f15bSPeter Brune     next->next      = ilink;
432eed5f15bSPeter Brune     ilink->previous = next;
433eed5f15bSPeter Brune   }
434eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
435eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
436eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
437eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
4388f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
439eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
4408f626970SPeter Brune   ilink->dmp = 1.0;
441*90a8ba9bSPeter Brune   jac->nsnes++;
442eed5f15bSPeter Brune   PetscFunctionReturn(0);
443eed5f15bSPeter Brune }
444eed5f15bSPeter Brune 
445eed5f15bSPeter Brune #undef __FUNCT__
446eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite"
447eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
448eed5f15bSPeter Brune {
449eed5f15bSPeter Brune   SNES_Composite     *jac;
450eed5f15bSPeter Brune   SNES_CompositeLink next;
451eed5f15bSPeter Brune   PetscInt         i;
452eed5f15bSPeter Brune 
453eed5f15bSPeter Brune   PetscFunctionBegin;
454eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
455eed5f15bSPeter Brune   next = jac->head;
456eed5f15bSPeter Brune   for (i=0; i<n; i++) {
457eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
458eed5f15bSPeter Brune     next = next->next;
459eed5f15bSPeter Brune   }
460eed5f15bSPeter Brune   *subsnes = next->snes;
461eed5f15bSPeter Brune   PetscFunctionReturn(0);
462eed5f15bSPeter Brune }
463eed5f15bSPeter Brune 
464eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
465eed5f15bSPeter Brune #undef __FUNCT__
466eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType"
467e31ab4f9SPeter Brune /*@C
468eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
469eed5f15bSPeter Brune 
470eed5f15bSPeter Brune    Logically Collective on SNES
471eed5f15bSPeter Brune 
472eed5f15bSPeter Brune    Input Parameter:
473eed5f15bSPeter Brune +  snes - the preconditioner context
474eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
475eed5f15bSPeter Brune 
476eed5f15bSPeter Brune    Options Database Key:
477eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
478eed5f15bSPeter Brune 
479eed5f15bSPeter Brune    Level: Developer
480eed5f15bSPeter Brune 
481eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
482eed5f15bSPeter Brune @*/
483eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
484eed5f15bSPeter Brune {
485eed5f15bSPeter Brune   PetscErrorCode ierr;
486eed5f15bSPeter Brune 
487eed5f15bSPeter Brune   PetscFunctionBegin;
488eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
489eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
490eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
491eed5f15bSPeter Brune   PetscFunctionReturn(0);
492eed5f15bSPeter Brune }
493eed5f15bSPeter Brune 
494eed5f15bSPeter Brune #undef __FUNCT__
495eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES"
496eed5f15bSPeter Brune /*@C
497eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
498eed5f15bSPeter Brune 
499eed5f15bSPeter Brune    Collective on SNES
500eed5f15bSPeter Brune 
501eed5f15bSPeter Brune    Input Parameters:
502eed5f15bSPeter Brune +  snes - the preconditioner context
503eed5f15bSPeter Brune -  type - the type of the new preconditioner
504eed5f15bSPeter Brune 
505eed5f15bSPeter Brune    Level: Developer
506eed5f15bSPeter Brune 
507eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add
508eed5f15bSPeter Brune @*/
509eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
510eed5f15bSPeter Brune {
511eed5f15bSPeter Brune   PetscErrorCode ierr;
512eed5f15bSPeter Brune 
513eed5f15bSPeter Brune   PetscFunctionBegin;
514eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
515eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
516eed5f15bSPeter Brune   PetscFunctionReturn(0);
517eed5f15bSPeter Brune }
518eed5f15bSPeter Brune #undef __FUNCT__
519eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES"
520eed5f15bSPeter Brune /*@
521eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
522eed5f15bSPeter Brune 
523eed5f15bSPeter Brune    Not Collective
524eed5f15bSPeter Brune 
525eed5f15bSPeter Brune    Input Parameter:
526eed5f15bSPeter Brune +  snes - the preconditioner context
527eed5f15bSPeter Brune -  n - the number of the snes requested
528eed5f15bSPeter Brune 
529eed5f15bSPeter Brune    Output Parameters:
530eed5f15bSPeter Brune .  subsnes - the SNES requested
531eed5f15bSPeter Brune 
532eed5f15bSPeter Brune    Level: Developer
533eed5f15bSPeter Brune 
534eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
535eed5f15bSPeter Brune 
536eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
537eed5f15bSPeter Brune @*/
538eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
539eed5f15bSPeter Brune {
540eed5f15bSPeter Brune   PetscErrorCode ierr;
541eed5f15bSPeter Brune 
542eed5f15bSPeter Brune   PetscFunctionBegin;
543eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
544eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
545eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
546eed5f15bSPeter Brune   PetscFunctionReturn(0);
547eed5f15bSPeter Brune }
548eed5f15bSPeter Brune 
549eed5f15bSPeter Brune #undef __FUNCT__
5508f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite"
5518f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
5528f626970SPeter Brune {
5538f626970SPeter Brune   SNES_Composite     *jac;
5548f626970SPeter Brune   SNES_CompositeLink next;
5558f626970SPeter Brune   PetscInt         i;
5568f626970SPeter Brune 
5578f626970SPeter Brune   PetscFunctionBegin;
5588f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
5598f626970SPeter Brune   next = jac->head;
5608f626970SPeter Brune   for (i=0; i<n; i++) {
5618f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
5628f626970SPeter Brune     next = next->next;
5638f626970SPeter Brune   }
5648f626970SPeter Brune   next->dmp = dmp;
5658f626970SPeter Brune   PetscFunctionReturn(0);
5668f626970SPeter Brune }
5678f626970SPeter Brune 
5688f626970SPeter Brune #undef __FUNCT__
5698f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping"
5708f626970SPeter Brune /*@
5718f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
5728f626970SPeter Brune 
5738f626970SPeter Brune    Not Collective
5748f626970SPeter Brune 
5758f626970SPeter Brune    Input Parameter:
5768f626970SPeter Brune +  snes - the preconditioner context
5778f626970SPeter Brune .  n - the number of the snes requested
5788f626970SPeter Brune -  dmp - the damping
5798f626970SPeter Brune 
5808f626970SPeter Brune    Level: Developer
5818f626970SPeter Brune 
5828f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
5838f626970SPeter Brune 
5848f626970SPeter Brune .seealso: SNESCompositeAddSNES()
5858f626970SPeter Brune @*/
5868f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
5878f626970SPeter Brune {
5888f626970SPeter Brune   PetscErrorCode ierr;
5898f626970SPeter Brune 
5908f626970SPeter Brune   PetscFunctionBegin;
5918f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5928f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
5938f626970SPeter Brune   PetscFunctionReturn(0);
5948f626970SPeter Brune }
5958f626970SPeter Brune 
5968f626970SPeter Brune #undef __FUNCT__
597eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite"
598eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes)
599eed5f15bSPeter Brune {
600eed5f15bSPeter Brune   Vec            F;
601eed5f15bSPeter Brune   Vec            X;
602eed5f15bSPeter Brune   Vec            B;
603eed5f15bSPeter Brune   PetscInt       i;
604eed5f15bSPeter Brune   PetscReal      fnorm = 0.0;
605eed5f15bSPeter Brune   PetscErrorCode ierr;
606eed5f15bSPeter Brune   SNESNormType   normtype;
607eed5f15bSPeter Brune   SNES_Composite *comp = (SNES_Composite*)snes->data;
608eed5f15bSPeter Brune 
609eed5f15bSPeter Brune   PetscFunctionBegin;
610eed5f15bSPeter Brune   X = snes->vec_sol;
611eed5f15bSPeter Brune   F = snes->vec_func;
612eed5f15bSPeter Brune   B = snes->vec_rhs;
613eed5f15bSPeter Brune 
614eed5f15bSPeter Brune   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
615eed5f15bSPeter Brune   snes->iter   = 0;
616eed5f15bSPeter Brune   snes->norm   = 0.;
617eed5f15bSPeter Brune   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
618eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
619eed5f15bSPeter Brune   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
620eed5f15bSPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
621eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
622eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
623eed5f15bSPeter Brune       if (snes->domainerror) {
624eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
625eed5f15bSPeter Brune         PetscFunctionReturn(0);
626eed5f15bSPeter Brune       }
627eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
628eed5f15bSPeter Brune 
629eed5f15bSPeter Brune     /* convergence test */
630eed5f15bSPeter Brune     if (!snes->norm_init_set) {
631eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
632eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
633eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
634eed5f15bSPeter Brune         PetscFunctionReturn(0);
635eed5f15bSPeter Brune       }
636eed5f15bSPeter Brune     } else {
637eed5f15bSPeter Brune       fnorm               = snes->norm_init;
638eed5f15bSPeter Brune       snes->norm_init_set = PETSC_FALSE;
639eed5f15bSPeter Brune     }
640eed5f15bSPeter Brune     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
641eed5f15bSPeter Brune     snes->iter = 0;
642eed5f15bSPeter Brune     snes->norm = fnorm;
643eed5f15bSPeter Brune     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
644eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
645eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
646eed5f15bSPeter Brune 
647eed5f15bSPeter Brune     /* set parameter for default relative tolerance convergence test */
648eed5f15bSPeter Brune     snes->ttol = fnorm*snes->rtol;
649eed5f15bSPeter Brune 
650eed5f15bSPeter Brune     /* test convergence */
651eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
652eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
653eed5f15bSPeter Brune   } else {
654eed5f15bSPeter Brune     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
655eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
656eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
657eed5f15bSPeter Brune   }
658eed5f15bSPeter Brune 
659eed5f15bSPeter Brune   /* Call general purpose update function */
660eed5f15bSPeter Brune   if (snes->ops->update) {
661eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
662eed5f15bSPeter Brune   }
663eed5f15bSPeter Brune 
664eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
665eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
666eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
667eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
668eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
669*90a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
670*90a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
671eed5f15bSPeter Brune     } else {
672*90a8ba9bSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
673eed5f15bSPeter Brune     }
674eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
675eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
676eed5f15bSPeter Brune       if (snes->domainerror) {
677eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
678eed5f15bSPeter Brune         break;
679eed5f15bSPeter Brune       }
680eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
681eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
682eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
683eed5f15bSPeter Brune         break;
684eed5f15bSPeter Brune       }
685eed5f15bSPeter Brune     }
686eed5f15bSPeter Brune     /* Monitor convergence */
687eed5f15bSPeter Brune     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
688eed5f15bSPeter Brune     snes->iter = i+1;
689eed5f15bSPeter Brune     snes->norm = fnorm;
690eed5f15bSPeter Brune     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
691eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
692eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
693eed5f15bSPeter Brune     /* Test for convergence */
694eed5f15bSPeter Brune     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
695eed5f15bSPeter Brune     if (snes->reason) break;
696eed5f15bSPeter Brune     /* Call general purpose update function */
697eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
698eed5f15bSPeter Brune   }
699eed5f15bSPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
700eed5f15bSPeter Brune     if (i == snes->max_its) {
701eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
702eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
703eed5f15bSPeter Brune     }
704eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
705eed5f15bSPeter Brune   PetscFunctionReturn(0);
706eed5f15bSPeter Brune }
707eed5f15bSPeter Brune 
708eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
709eed5f15bSPeter Brune 
710eed5f15bSPeter Brune /*MC
711eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
712eed5f15bSPeter Brune 
713eed5f15bSPeter Brune    Options Database Keys:
714eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
715eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
716eed5f15bSPeter Brune 
717eed5f15bSPeter Brune    Level: intermediate
718eed5f15bSPeter Brune 
719eed5f15bSPeter Brune    Concepts: composing solvers
720eed5f15bSPeter Brune 
721eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
722eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
723eed5f15bSPeter Brune            SNESCompositeGetSNES()
724eed5f15bSPeter Brune 
725eed5f15bSPeter Brune M*/
726eed5f15bSPeter Brune 
727eed5f15bSPeter Brune #undef __FUNCT__
728eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite"
729eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
730eed5f15bSPeter Brune {
731eed5f15bSPeter Brune   PetscErrorCode ierr;
732eed5f15bSPeter Brune   SNES_Composite   *jac;
733eed5f15bSPeter Brune 
734eed5f15bSPeter Brune   PetscFunctionBegin;
735eed5f15bSPeter Brune   ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr);
736eed5f15bSPeter Brune 
737eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
738eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
739eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
740eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
741eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
742eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
743eed5f15bSPeter Brune 
744eed5f15bSPeter Brune   snes->data = (void*)jac;
745*90a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
746*90a8ba9bSPeter Brune   jac->Fes   = NULL;
747*90a8ba9bSPeter Brune   jac->Xes   = NULL;
748*90a8ba9bSPeter Brune   jac->nsnes = 0;
749eed5f15bSPeter Brune   jac->head  = 0;
750eed5f15bSPeter Brune 
751*90a8ba9bSPeter Brune   jac->h     = NULL;
752*90a8ba9bSPeter Brune   jac->s     = NULL;
753*90a8ba9bSPeter Brune   jac->beta  = NULL;
754*90a8ba9bSPeter Brune   jac->work  = NULL;
755*90a8ba9bSPeter Brune   jac->rwork = NULL;
756*90a8ba9bSPeter Brune 
7578f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
7588f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
7598f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
7608f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
761eed5f15bSPeter Brune   PetscFunctionReturn(0);
762eed5f15bSPeter Brune }
763eed5f15bSPeter Brune 
764