xref: /petsc/src/snes/impls/composite/snescomposite.c (revision 63cdc2bd36998bbfb989661ff0ef5402b58727af)
1eed5f15bSPeter Brune 
2eed5f15bSPeter Brune /*
3eed5f15bSPeter Brune       Defines a SNES that can consist of a collection of SNESes
4eed5f15bSPeter Brune */
5af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
690a8ba9bSPeter Brune #include <petscblaslapack.h>
7eed5f15bSPeter Brune 
890a8ba9bSPeter 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;
1490a8ba9bSPeter 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;
2190a8ba9bSPeter Brune   PetscInt           nsnes;
22eed5f15bSPeter Brune   SNESCompositeType  type;
23eed5f15bSPeter Brune   Vec                Xorig;
2490a8ba9bSPeter Brune 
2590a8ba9bSPeter Brune   /* context for ADDITIVEOPTIMAL */
2690a8ba9bSPeter Brune   Vec                *Xes,*Fes;      /* solution and residual vectors for the subsolvers */
279c2f3473SPeter Brune   PetscReal          *fnorms;        /* norms of the solutions */
2890a8ba9bSPeter Brune   PetscScalar        *h;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
299c2f3473SPeter Brune   PetscScalar        *g;             /* the dotproducts of the previous function with the candidate functions */
309c2f3473SPeter Brune   PetscBLASInt       n;              /* matrix dimension -- nsnes */
3190a8ba9bSPeter Brune   PetscBLASInt       nrhs;           /* the number of right hand sides */
3290a8ba9bSPeter Brune   PetscBLASInt       lda;            /* the padded matrix dimension */
3390a8ba9bSPeter Brune   PetscBLASInt       ldb;            /* the padded vector dimension */
3490a8ba9bSPeter Brune   PetscReal          *s;             /* the singular values */
359c2f3473SPeter Brune   PetscScalar        *beta;          /* the RHS and combination */
3690a8ba9bSPeter Brune   PetscReal          rcond;          /* the exit condition */
3790a8ba9bSPeter Brune   PetscBLASInt       rank;           /* the effective rank */
3890a8ba9bSPeter Brune   PetscScalar        *work;          /* the work vector */
3990a8ba9bSPeter Brune   PetscReal          *rwork;         /* the real work vector used for complex */
4090a8ba9bSPeter Brune   PetscBLASInt       lwork;          /* the size of the work vector */
4190a8ba9bSPeter Brune   PetscBLASInt       info;           /* the output condition */
4290a8ba9bSPeter Brune 
435e806d2eSPeter Brune   PetscReal          rtol;           /* restart tolerance for accepting the combination */
445e806d2eSPeter Brune   PetscReal          stol;           /* restart tolerance for the combination */
45eed5f15bSPeter Brune } SNES_Composite;
46eed5f15bSPeter Brune 
47eed5f15bSPeter Brune #undef __FUNCT__
48eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative"
49eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
50eed5f15bSPeter Brune {
51eed5f15bSPeter Brune   PetscErrorCode      ierr;
52eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
53eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
54eed5f15bSPeter Brune   Vec                 FSub;
55d5a53f06SPeter Brune   SNESConvergedReason reason;
56eed5f15bSPeter Brune 
57eed5f15bSPeter Brune   PetscFunctionBegin;
58eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
5972edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
60eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
61eed5f15bSPeter Brune   }
62eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
63d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
64d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
65d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
66d5a53f06SPeter Brune     PetscFunctionReturn(0);
67d5a53f06SPeter Brune   }
68eed5f15bSPeter Brune 
69eed5f15bSPeter Brune   while (next->next) {
70eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
7172edecb9SPeter Brune     if (next->snes->pcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
72eed5f15bSPeter Brune       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
73eed5f15bSPeter Brune       next = next->next;
74eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
75eed5f15bSPeter Brune     } else {
76eed5f15bSPeter Brune       next = next->next;
77eed5f15bSPeter Brune     }
78eed5f15bSPeter Brune     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
79d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
80d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
81d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
82d5a53f06SPeter Brune       PetscFunctionReturn(0);
83d5a53f06SPeter Brune     }
84eed5f15bSPeter Brune   }
85eed5f15bSPeter Brune   if (next->snes->pcside == PC_RIGHT) {
86eed5f15bSPeter Brune     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
87eed5f15bSPeter Brune     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
88d5a53f06SPeter Brune     if (fnorm) {
89*63cdc2bdSPatrick Farrell       if (snes->xl && snes->xu) {
90*63cdc2bdSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
91*63cdc2bdSPatrick Farrell       } else {
9271dbe336SPeter Brune         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
93*63cdc2bdSPatrick Farrell       }
94*63cdc2bdSPatrick Farrell 
95d5a53f06SPeter Brune       if (PetscIsInfOrNanReal(*fnorm)) {
96d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
97d5a53f06SPeter Brune         PetscFunctionReturn(0);
98d5a53f06SPeter Brune       }
99d5a53f06SPeter Brune     }
10072edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
101be95d8f1SBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
102d5a53f06SPeter Brune     if (snes->domainerror) {
103d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
104d5a53f06SPeter Brune       PetscFunctionReturn(0);
105d5a53f06SPeter Brune     }
106d5a53f06SPeter Brune     if (fnorm) {
107a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
108a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
109a6da83ecSPatrick Farrell       } else {
110d5a53f06SPeter Brune         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
111a6da83ecSPatrick Farrell       }
112a6da83ecSPatrick Farrell 
113d5a53f06SPeter Brune       if (PetscIsInfOrNanReal(*fnorm)) {
114d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
115d5a53f06SPeter Brune         PetscFunctionReturn(0);
116d5a53f06SPeter Brune       }
117d5a53f06SPeter Brune     }
118eed5f15bSPeter Brune   }
119eed5f15bSPeter Brune   PetscFunctionReturn(0);
120eed5f15bSPeter Brune }
121eed5f15bSPeter Brune 
122eed5f15bSPeter Brune #undef __FUNCT__
123eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive"
124eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
125eed5f15bSPeter Brune {
126eed5f15bSPeter Brune   PetscErrorCode      ierr;
127eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
128eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
129eed5f15bSPeter Brune   Vec                 Y,Xorig;
130d5a53f06SPeter Brune   SNESConvergedReason reason;
131eed5f15bSPeter Brune 
132eed5f15bSPeter Brune   PetscFunctionBegin;
133eed5f15bSPeter Brune   Y = snes->vec_sol_update;
134eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
135eed5f15bSPeter Brune   Xorig = jac->Xorig;
136302440fdSBarry Smith   ierr = VecCopy(X,Xorig);CHKERRQ(ierr);
137eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
13872edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
139eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
140eed5f15bSPeter Brune     while (next->next) {
141eed5f15bSPeter Brune       next = next->next;
142eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
143eed5f15bSPeter Brune     }
144eed5f15bSPeter Brune   }
145eed5f15bSPeter Brune   next = jac->head;
1468f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
147eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
148d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
149d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
150d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
151d5a53f06SPeter Brune     PetscFunctionReturn(0);
152d5a53f06SPeter Brune   }
153eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1548f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1558f626970SPeter Brune   while (next->next) {
1568f626970SPeter Brune     next = next->next;
1578f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1588f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
159d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
160d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
161d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
162d5a53f06SPeter Brune       PetscFunctionReturn(0);
163d5a53f06SPeter Brune     }
1648f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1658f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
166eed5f15bSPeter Brune   }
16772edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
168eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
169a6da83ecSPatrick Farrell     if (fnorm) {
170a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
171a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
172a6da83ecSPatrick Farrell       } else {
173a6da83ecSPatrick Farrell         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
174a6da83ecSPatrick Farrell       }
175a6da83ecSPatrick Farrell 
176a6da83ecSPatrick Farrell       if (PetscIsInfOrNanReal(*fnorm)) {
177a6da83ecSPatrick Farrell         snes->reason = SNES_DIVERGED_FNORM_NAN;
178a6da83ecSPatrick Farrell         PetscFunctionReturn(0);
179a6da83ecSPatrick Farrell       }
180a6da83ecSPatrick Farrell     }
181eed5f15bSPeter Brune   }
182eed5f15bSPeter Brune   PetscFunctionReturn(0);
183eed5f15bSPeter Brune }
18490a8ba9bSPeter Brune 
18590a8ba9bSPeter Brune #undef __FUNCT__
18690a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal"
18790a8ba9bSPeter Brune /* approximately solve the overdetermined system:
18890a8ba9bSPeter Brune 
18990a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
19090a8ba9bSPeter Brune  \alpha_i                      = 1
19190a8ba9bSPeter Brune 
19290a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
19390a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
19490a8ba9bSPeter Brune 
19590a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
19690a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
19790a8ba9bSPeter Brune  */
19890a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
19990a8ba9bSPeter Brune {
20090a8ba9bSPeter Brune   PetscErrorCode      ierr;
20190a8ba9bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
20290a8ba9bSPeter Brune   SNES_CompositeLink  next = jac->head;
20390a8ba9bSPeter Brune   Vec                 *Xes = jac->Xes,*Fes = jac->Fes;
20490a8ba9bSPeter Brune   PetscInt            i,j;
2055e806d2eSPeter Brune   PetscScalar         tot,total,ftf;
2065e806d2eSPeter Brune   PetscReal           min_fnorm;
2075e806d2eSPeter Brune   PetscInt            min_i;
208d5a53f06SPeter Brune   SNESConvergedReason reason;
20990a8ba9bSPeter Brune 
21090a8ba9bSPeter Brune   PetscFunctionBegin;
21190a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
21258bdfa14SPeter Brune 
21358bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
21458bdfa14SPeter Brune     next = jac->head;
21558bdfa14SPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
21658bdfa14SPeter Brune     while (next->next) {
21758bdfa14SPeter Brune       next = next->next;
21858bdfa14SPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
21958bdfa14SPeter Brune     }
22058bdfa14SPeter Brune   }
22158bdfa14SPeter Brune 
22290a8ba9bSPeter Brune   next = jac->head;
22390a8ba9bSPeter Brune   i = 0;
22490a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
22590a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
226d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
227d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
228d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
229d5a53f06SPeter Brune     PetscFunctionReturn(0);
230d5a53f06SPeter Brune   }
23190a8ba9bSPeter Brune   while (next->next) {
23290a8ba9bSPeter Brune     i++;
23390a8ba9bSPeter Brune     next = next->next;
23490a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
23590a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
236d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
237d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
238d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
239d5a53f06SPeter Brune       PetscFunctionReturn(0);
240d5a53f06SPeter Brune     }
24190a8ba9bSPeter Brune   }
24290a8ba9bSPeter Brune 
24390a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
24490a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
24590a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2469c2f3473SPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
24790a8ba9bSPeter Brune     }
2489c2f3473SPeter Brune     ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
24990a8ba9bSPeter Brune   }
2505e806d2eSPeter Brune 
25190a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
25290a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2539c2f3473SPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
2549c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
25590a8ba9bSPeter Brune     }
2569c2f3473SPeter Brune     ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
25790a8ba9bSPeter Brune   }
25890a8ba9bSPeter Brune 
2599c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
26090a8ba9bSPeter Brune 
26190a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
26290a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
2639c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
26490a8ba9bSPeter Brune     }
26590a8ba9bSPeter Brune   }
26690a8ba9bSPeter Brune 
26790a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
2689c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
2699c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
27090a8ba9bSPeter Brune     }
2719c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2729c2f3473SPeter Brune   }
27390a8ba9bSPeter Brune 
27490a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
27590a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
27690a8ba9bSPeter Brune #else
27790a8ba9bSPeter Brune   jac->info  = 0;
27890a8ba9bSPeter Brune   jac->rcond = -1.;
27990a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
28090a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
2819c2f3473SPeter Brune   PetscStackCall("LAPACKgelss",LAPACKgelss_(&jac->n,&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));
28290a8ba9bSPeter Brune #else
2839c2f3473SPeter Brune   PetscStackCall("LAPACKgelss",LAPACKgelss_(&jac->n,&jac->n,&jac->nrhs,jac->h,&jac->lda,jac->beta,&jac->lda,jac->s,&jac->rcond,&jac->rank,jac->work,&jac->lwork,&jac->info));
28490a8ba9bSPeter Brune #endif
28590a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
28690a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
28790a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
28890a8ba9bSPeter Brune #endif
2899c2f3473SPeter Brune   tot = 0.;
2905e806d2eSPeter Brune   total = 0.;
29190a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
29290a8ba9bSPeter Brune     if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
293c783eb9eSBarry Smith     ierr = PetscInfo2(snes,"%D: %g\n",i,(double)PetscRealPart(jac->beta[i]));CHKERRQ(ierr);
2949c2f3473SPeter Brune     tot += jac->beta[i];
2955e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
29690a8ba9bSPeter Brune   }
2979c2f3473SPeter Brune   ierr = VecScale(X,(1. - tot));CHKERRQ(ierr);
29890a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
29990a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
300a6da83ecSPatrick Farrell 
301a6da83ecSPatrick Farrell   if (snes->xl && snes->xu) {
302a6da83ecSPatrick Farrell     ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
303a6da83ecSPatrick Farrell   } else {
3049c2f3473SPeter Brune     ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
305a6da83ecSPatrick Farrell   }
30690a8ba9bSPeter Brune 
3075e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
3085e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
3095e806d2eSPeter Brune   min_i     = 0;
3109c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
3115e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
3125e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
3135e806d2eSPeter Brune       min_i     = i;
3149c2f3473SPeter Brune     }
3159c2f3473SPeter Brune   }
3165e806d2eSPeter Brune 
3175e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
3185e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
31958bdfa14SPeter Brune     ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr);
32058bdfa14SPeter Brune     ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr);
3215e806d2eSPeter Brune     *fnorm = min_fnorm;
3225e806d2eSPeter Brune   }
32390a8ba9bSPeter Brune   PetscFunctionReturn(0);
32490a8ba9bSPeter Brune }
32590a8ba9bSPeter Brune 
326eed5f15bSPeter Brune #undef __FUNCT__
327eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite"
328eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
329eed5f15bSPeter Brune {
330eed5f15bSPeter Brune   PetscErrorCode     ierr;
331dd515a93SPeter Brune   DM                 dm;
332eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
333eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
33490a8ba9bSPeter Brune   PetscInt           n=0,i;
33590a8ba9bSPeter Brune   Vec                F;
336eed5f15bSPeter Brune 
337eed5f15bSPeter Brune   PetscFunctionBegin;
338dd515a93SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
339*63cdc2bdSPatrick Farrell 
340*63cdc2bdSPatrick Farrell   if (snes->ops->computevariablebounds) {
341*63cdc2bdSPatrick Farrell     /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */
342*63cdc2bdSPatrick Farrell     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
343*63cdc2bdSPatrick Farrell     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
344*63cdc2bdSPatrick Farrell     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
345*63cdc2bdSPatrick Farrell   }
346*63cdc2bdSPatrick Farrell 
347eed5f15bSPeter Brune   while (next) {
34890a8ba9bSPeter Brune     n++;
349dd515a93SPeter Brune     ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr);
350eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
351*63cdc2bdSPatrick Farrell     ierr = SNESSetApplicationContext(next->snes, snes->user);CHKERRQ(ierr);
352*63cdc2bdSPatrick Farrell     if (snes->xl && snes->xu) {
353*63cdc2bdSPatrick Farrell       if (snes->ops->computevariablebounds) {
354*63cdc2bdSPatrick Farrell         ierr = SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds);CHKERRQ(ierr);
355*63cdc2bdSPatrick Farrell       } else {
356*63cdc2bdSPatrick Farrell         ierr = SNESVISetVariableBounds(next->snes,snes->xl,snes->xu);CHKERRQ(ierr);
357*63cdc2bdSPatrick Farrell       }
358*63cdc2bdSPatrick Farrell     }
359*63cdc2bdSPatrick Farrell 
360eed5f15bSPeter Brune     next = next->next;
361eed5f15bSPeter Brune   }
36290a8ba9bSPeter Brune   jac->nsnes = n;
36390a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
36490a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
36590a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
366854ce69bSBarry Smith     ierr = PetscMalloc1(n,&jac->Fes);CHKERRQ(ierr);
367854ce69bSBarry Smith     ierr = PetscMalloc1(n,&jac->fnorms);CHKERRQ(ierr);
36890a8ba9bSPeter Brune     next = jac->head;
36990a8ba9bSPeter Brune     i = 0;
37090a8ba9bSPeter Brune     while (next) {
37190a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
37290a8ba9bSPeter Brune       jac->Fes[i] = F;
37390a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
37490a8ba9bSPeter Brune       next = next->next;
37590a8ba9bSPeter Brune       i++;
37690a8ba9bSPeter Brune     }
37790a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
37890a8ba9bSPeter Brune     jac->nrhs  = 1;
3799c2f3473SPeter Brune     jac->lda   = jac->nsnes;
38090a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
38190a8ba9bSPeter Brune     jac->n     = jac->nsnes;
38290a8ba9bSPeter Brune 
383785e854fSJed Brown     ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr);
384785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr);
385785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr);
386785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr);
3879c2f3473SPeter Brune     jac->lwork = 12*jac->n;
38890a8ba9bSPeter Brune #if PETSC_USE_COMPLEX
389854ce69bSBarry Smith     ierr = PetscMalloc1(jac->lwork,&jac->rwork);CHKERRQ(ierr);
39090a8ba9bSPeter Brune #endif
391854ce69bSBarry Smith     ierr = PetscMalloc1(jac->lwork,&jac->work);CHKERRQ(ierr);
39290a8ba9bSPeter Brune   }
39390a8ba9bSPeter Brune 
394eed5f15bSPeter Brune   PetscFunctionReturn(0);
395eed5f15bSPeter Brune }
396eed5f15bSPeter Brune 
397eed5f15bSPeter Brune #undef __FUNCT__
398eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite"
399eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
400eed5f15bSPeter Brune {
401eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
402eed5f15bSPeter Brune   PetscErrorCode   ierr;
403eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
404eed5f15bSPeter Brune 
405eed5f15bSPeter Brune   PetscFunctionBegin;
406eed5f15bSPeter Brune   while (next) {
407eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
408eed5f15bSPeter Brune     next = next->next;
409eed5f15bSPeter Brune   }
410eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
41190a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
41290a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
4139c2f3473SPeter Brune   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
41490a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
41590a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
4169c2f3473SPeter Brune   ierr = PetscFree(jac->g);CHKERRQ(ierr);
41790a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
41890a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
41990a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
420eed5f15bSPeter Brune   PetscFunctionReturn(0);
421eed5f15bSPeter Brune }
422eed5f15bSPeter Brune 
423eed5f15bSPeter Brune #undef __FUNCT__
424eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite"
425eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
426eed5f15bSPeter Brune {
427eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
428eed5f15bSPeter Brune   PetscErrorCode   ierr;
429eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
430eed5f15bSPeter Brune 
431eed5f15bSPeter Brune   PetscFunctionBegin;
432eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
433eed5f15bSPeter Brune   while (next) {
434eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
435eed5f15bSPeter Brune     next_tmp = next;
436eed5f15bSPeter Brune     next     = next->next;
437eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
438eed5f15bSPeter Brune   }
439eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
440eed5f15bSPeter Brune   PetscFunctionReturn(0);
441eed5f15bSPeter Brune }
442eed5f15bSPeter Brune 
443eed5f15bSPeter Brune #undef __FUNCT__
444eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite"
4458c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(PetscOptions *PetscOptionsObject,SNES snes)
446eed5f15bSPeter Brune {
447eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
448eed5f15bSPeter Brune   PetscErrorCode     ierr;
449eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
450eed5f15bSPeter Brune   SNES_CompositeLink next;
451eed5f15bSPeter Brune   char               *sneses[8];
4528f626970SPeter Brune   PetscReal          dmps[8];
453eed5f15bSPeter Brune   PetscBool          flg;
454eed5f15bSPeter Brune 
455eed5f15bSPeter Brune   PetscFunctionBegin;
456e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");CHKERRQ(ierr);
457eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
458eed5f15bSPeter Brune   if (flg) {
459eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
460eed5f15bSPeter Brune   }
461eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
462eed5f15bSPeter Brune   if (flg) {
463eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
464eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
465eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
466eed5f15bSPeter Brune     }
467eed5f15bSPeter Brune   }
4688f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
4698f626970SPeter Brune   if (flg) {
4708f626970SPeter Brune     for (i=0; i<nmax; i++) {
4718f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
4728f626970SPeter Brune     }
4738f626970SPeter Brune   }
4745e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr);
4755e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr);
476eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
477eed5f15bSPeter Brune 
478eed5f15bSPeter Brune   next = jac->head;
479eed5f15bSPeter Brune   while (next) {
480eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
481eed5f15bSPeter Brune     next = next->next;
482eed5f15bSPeter Brune   }
483eed5f15bSPeter Brune   PetscFunctionReturn(0);
484eed5f15bSPeter Brune }
485eed5f15bSPeter Brune 
486eed5f15bSPeter Brune #undef __FUNCT__
487eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite"
488eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
489eed5f15bSPeter Brune {
490eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
491eed5f15bSPeter Brune   PetscErrorCode   ierr;
492eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
493eed5f15bSPeter Brune   PetscBool        iascii;
494eed5f15bSPeter Brune 
495eed5f15bSPeter Brune   PetscFunctionBegin;
496eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
497eed5f15bSPeter Brune   if (iascii) {
498eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
499eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
500eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
501eed5f15bSPeter Brune   }
502eed5f15bSPeter Brune   if (iascii) {
503eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
504eed5f15bSPeter Brune   }
505eed5f15bSPeter Brune   while (next) {
506eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
507eed5f15bSPeter Brune     next = next->next;
508eed5f15bSPeter Brune   }
509eed5f15bSPeter Brune   if (iascii) {
510eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
511eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
512eed5f15bSPeter Brune   }
513eed5f15bSPeter Brune   PetscFunctionReturn(0);
514eed5f15bSPeter Brune }
515eed5f15bSPeter Brune 
516eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
517eed5f15bSPeter Brune 
518eed5f15bSPeter Brune #undef __FUNCT__
519eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite"
520eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
521eed5f15bSPeter Brune {
522eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
523eed5f15bSPeter Brune 
524eed5f15bSPeter Brune   PetscFunctionBegin;
525eed5f15bSPeter Brune   jac->type = type;
526eed5f15bSPeter Brune   PetscFunctionReturn(0);
527eed5f15bSPeter Brune }
528eed5f15bSPeter Brune 
529eed5f15bSPeter Brune #undef __FUNCT__
530eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite"
531eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
532eed5f15bSPeter Brune {
533eed5f15bSPeter Brune   SNES_Composite     *jac;
534eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
535eed5f15bSPeter Brune   PetscErrorCode   ierr;
536eed5f15bSPeter Brune   PetscInt         cnt = 0;
537eed5f15bSPeter Brune   const char       *prefix;
538eed5f15bSPeter Brune   char             newprefix[8];
539eed5f15bSPeter Brune   DM               dm;
540eed5f15bSPeter Brune 
541eed5f15bSPeter Brune   PetscFunctionBegin;
542b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ilink);CHKERRQ(ierr);
543eed5f15bSPeter Brune   ilink->next = 0;
544eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
545eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
546eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
547eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
548cf5b3eb5SPeter Brune   ierr        = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr);
549eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
550eed5f15bSPeter Brune   next = jac->head;
551eed5f15bSPeter Brune   if (!next) {
552eed5f15bSPeter Brune     jac->head       = ilink;
553eed5f15bSPeter Brune     ilink->previous = NULL;
554eed5f15bSPeter Brune   } else {
555eed5f15bSPeter Brune     cnt++;
556eed5f15bSPeter Brune     while (next->next) {
557eed5f15bSPeter Brune       next = next->next;
558eed5f15bSPeter Brune       cnt++;
559eed5f15bSPeter Brune     }
560eed5f15bSPeter Brune     next->next      = ilink;
561eed5f15bSPeter Brune     ilink->previous = next;
562eed5f15bSPeter Brune   }
563eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
564eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
565eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
566eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
5678f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
568eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
56972edecb9SPeter Brune   ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
570*63cdc2bdSPatrick Farrell 
5718f626970SPeter Brune   ilink->dmp = 1.0;
57290a8ba9bSPeter Brune   jac->nsnes++;
573eed5f15bSPeter Brune   PetscFunctionReturn(0);
574eed5f15bSPeter Brune }
575eed5f15bSPeter Brune 
576eed5f15bSPeter Brune #undef __FUNCT__
577eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite"
578eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
579eed5f15bSPeter Brune {
580eed5f15bSPeter Brune   SNES_Composite     *jac;
581eed5f15bSPeter Brune   SNES_CompositeLink next;
582eed5f15bSPeter Brune   PetscInt         i;
583eed5f15bSPeter Brune 
584eed5f15bSPeter Brune   PetscFunctionBegin;
585eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
586eed5f15bSPeter Brune   next = jac->head;
587eed5f15bSPeter Brune   for (i=0; i<n; i++) {
588eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
589eed5f15bSPeter Brune     next = next->next;
590eed5f15bSPeter Brune   }
591eed5f15bSPeter Brune   *subsnes = next->snes;
592eed5f15bSPeter Brune   PetscFunctionReturn(0);
593eed5f15bSPeter Brune }
594eed5f15bSPeter Brune 
595eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
596eed5f15bSPeter Brune #undef __FUNCT__
597eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType"
598e31ab4f9SPeter Brune /*@C
599eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
600eed5f15bSPeter Brune 
601eed5f15bSPeter Brune    Logically Collective on SNES
602eed5f15bSPeter Brune 
603eed5f15bSPeter Brune    Input Parameter:
604eed5f15bSPeter Brune +  snes - the preconditioner context
605eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
606eed5f15bSPeter Brune 
607eed5f15bSPeter Brune    Options Database Key:
608eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
609eed5f15bSPeter Brune 
610eed5f15bSPeter Brune    Level: Developer
611eed5f15bSPeter Brune 
612eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
613eed5f15bSPeter Brune @*/
614eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
615eed5f15bSPeter Brune {
616eed5f15bSPeter Brune   PetscErrorCode ierr;
617eed5f15bSPeter Brune 
618eed5f15bSPeter Brune   PetscFunctionBegin;
619eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
620eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
621eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
622eed5f15bSPeter Brune   PetscFunctionReturn(0);
623eed5f15bSPeter Brune }
624eed5f15bSPeter Brune 
625eed5f15bSPeter Brune #undef __FUNCT__
626eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES"
627eed5f15bSPeter Brune /*@C
628eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
629eed5f15bSPeter Brune 
630eed5f15bSPeter Brune    Collective on SNES
631eed5f15bSPeter Brune 
632eed5f15bSPeter Brune    Input Parameters:
633eed5f15bSPeter Brune +  snes - the preconditioner context
634eed5f15bSPeter Brune -  type - the type of the new preconditioner
635eed5f15bSPeter Brune 
636eed5f15bSPeter Brune    Level: Developer
637eed5f15bSPeter Brune 
638eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add
639eed5f15bSPeter Brune @*/
640eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
641eed5f15bSPeter Brune {
642eed5f15bSPeter Brune   PetscErrorCode ierr;
643eed5f15bSPeter Brune 
644eed5f15bSPeter Brune   PetscFunctionBegin;
645eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
646eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
647eed5f15bSPeter Brune   PetscFunctionReturn(0);
648eed5f15bSPeter Brune }
649eed5f15bSPeter Brune #undef __FUNCT__
650eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES"
651eed5f15bSPeter Brune /*@
652eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
653eed5f15bSPeter Brune 
654eed5f15bSPeter Brune    Not Collective
655eed5f15bSPeter Brune 
656eed5f15bSPeter Brune    Input Parameter:
657eed5f15bSPeter Brune +  snes - the preconditioner context
658eed5f15bSPeter Brune -  n - the number of the snes requested
659eed5f15bSPeter Brune 
660eed5f15bSPeter Brune    Output Parameters:
661eed5f15bSPeter Brune .  subsnes - the SNES requested
662eed5f15bSPeter Brune 
663eed5f15bSPeter Brune    Level: Developer
664eed5f15bSPeter Brune 
665eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
666eed5f15bSPeter Brune 
667eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
668eed5f15bSPeter Brune @*/
669eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
670eed5f15bSPeter Brune {
671eed5f15bSPeter Brune   PetscErrorCode ierr;
672eed5f15bSPeter Brune 
673eed5f15bSPeter Brune   PetscFunctionBegin;
674eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
675eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
676eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
677eed5f15bSPeter Brune   PetscFunctionReturn(0);
678eed5f15bSPeter Brune }
679eed5f15bSPeter Brune 
680eed5f15bSPeter Brune #undef __FUNCT__
6816b2b7f7bSPatrick Farrell #define __FUNCT__ "SNESCompositeGetNumber"
6826b2b7f7bSPatrick Farrell /*@
6836b2b7f7bSPatrick Farrell    SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES.
6846b2b7f7bSPatrick Farrell 
6856b2b7f7bSPatrick Farrell    Logically Collective on SNES
6866b2b7f7bSPatrick Farrell 
6876b2b7f7bSPatrick Farrell    Input Parameter:
6886b2b7f7bSPatrick Farrell    snes - the preconditioner context
6896b2b7f7bSPatrick Farrell 
6906b2b7f7bSPatrick Farrell    Output Parameter:
6916b2b7f7bSPatrick Farrell    n - the number of subsolvers
6926b2b7f7bSPatrick Farrell 
6936b2b7f7bSPatrick Farrell    Level: Developer
6946b2b7f7bSPatrick Farrell 
6956b2b7f7bSPatrick Farrell .keywords: SNES, composite preconditioner
6966b2b7f7bSPatrick Farrell @*/
6976b2b7f7bSPatrick Farrell PetscErrorCode  SNESCompositeGetNumber(SNES snes,PetscInt *n)
6986b2b7f7bSPatrick Farrell {
6996b2b7f7bSPatrick Farrell   SNES_Composite     *jac;
7006b2b7f7bSPatrick Farrell   SNES_CompositeLink next;
7016b2b7f7bSPatrick Farrell 
7026b2b7f7bSPatrick Farrell   PetscFunctionBegin;
7036b2b7f7bSPatrick Farrell   jac  = (SNES_Composite*)snes->data;
7046b2b7f7bSPatrick Farrell   next = jac->head;
7056b2b7f7bSPatrick Farrell 
7066b2b7f7bSPatrick Farrell   *n = 0;
7076b2b7f7bSPatrick Farrell   while (next) {
7086b2b7f7bSPatrick Farrell     *n = *n + 1;
7096b2b7f7bSPatrick Farrell     next = next->next;
7106b2b7f7bSPatrick Farrell   }
7116b2b7f7bSPatrick Farrell   PetscFunctionReturn(0);
7126b2b7f7bSPatrick Farrell }
7136b2b7f7bSPatrick Farrell 
7146b2b7f7bSPatrick Farrell #undef __FUNCT__
7158f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite"
7168f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
7178f626970SPeter Brune {
7188f626970SPeter Brune   SNES_Composite     *jac;
7198f626970SPeter Brune   SNES_CompositeLink next;
7208f626970SPeter Brune   PetscInt         i;
7218f626970SPeter Brune 
7228f626970SPeter Brune   PetscFunctionBegin;
7238f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
7248f626970SPeter Brune   next = jac->head;
7258f626970SPeter Brune   for (i=0; i<n; i++) {
7268f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
7278f626970SPeter Brune     next = next->next;
7288f626970SPeter Brune   }
7298f626970SPeter Brune   next->dmp = dmp;
7308f626970SPeter Brune   PetscFunctionReturn(0);
7318f626970SPeter Brune }
7328f626970SPeter Brune 
7338f626970SPeter Brune #undef __FUNCT__
7348f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping"
7358f626970SPeter Brune /*@
7368f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
7378f626970SPeter Brune 
7388f626970SPeter Brune    Not Collective
7398f626970SPeter Brune 
7408f626970SPeter Brune    Input Parameter:
7418f626970SPeter Brune +  snes - the preconditioner context
7428f626970SPeter Brune .  n - the number of the snes requested
7438f626970SPeter Brune -  dmp - the damping
7448f626970SPeter Brune 
7458f626970SPeter Brune    Level: Developer
7468f626970SPeter Brune 
7478f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
7488f626970SPeter Brune 
7498f626970SPeter Brune .seealso: SNESCompositeAddSNES()
7508f626970SPeter Brune @*/
7518f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
7528f626970SPeter Brune {
7538f626970SPeter Brune   PetscErrorCode ierr;
7548f626970SPeter Brune 
7558f626970SPeter Brune   PetscFunctionBegin;
7568f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7578f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
7588f626970SPeter Brune   PetscFunctionReturn(0);
7598f626970SPeter Brune }
7608f626970SPeter Brune 
7618f626970SPeter Brune #undef __FUNCT__
762eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite"
763eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes)
764eed5f15bSPeter Brune {
765eed5f15bSPeter Brune   Vec            F;
766eed5f15bSPeter Brune   Vec            X;
767eed5f15bSPeter Brune   Vec            B;
768eed5f15bSPeter Brune   PetscInt       i;
769b22975d2SPatrick Farrell   PetscReal      fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
770eed5f15bSPeter Brune   PetscErrorCode ierr;
77172edecb9SPeter Brune   SNESNormSchedule normtype;
772eed5f15bSPeter Brune   SNES_Composite *comp = (SNES_Composite*)snes->data;
773eed5f15bSPeter Brune 
774eed5f15bSPeter Brune   PetscFunctionBegin;
775eed5f15bSPeter Brune   X = snes->vec_sol;
776eed5f15bSPeter Brune   F = snes->vec_func;
777eed5f15bSPeter Brune   B = snes->vec_rhs;
778eed5f15bSPeter Brune 
779e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
780eed5f15bSPeter Brune   snes->iter   = 0;
781eed5f15bSPeter Brune   snes->norm   = 0.;
782e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
783b22975d2SPatrick Farrell   ierr         = SNESSetWorkVecs(snes, 1);CHKERRQ(ierr);
784eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
78572edecb9SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr);
78672edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
787eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
788eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
789eed5f15bSPeter Brune       if (snes->domainerror) {
790eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
791eed5f15bSPeter Brune         PetscFunctionReturn(0);
792eed5f15bSPeter Brune       }
793eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
794eed5f15bSPeter Brune 
795a6da83ecSPatrick Farrell     if (snes->xl && snes->xu) {
796a6da83ecSPatrick Farrell       ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr);
797a6da83ecSPatrick Farrell     } else {
798eed5f15bSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
799a6da83ecSPatrick Farrell     }
800a6da83ecSPatrick Farrell 
801eed5f15bSPeter Brune     if (PetscIsInfOrNanReal(fnorm)) {
802eed5f15bSPeter Brune       snes->reason = SNES_DIVERGED_FNORM_NAN;
803eed5f15bSPeter Brune       PetscFunctionReturn(0);
804eed5f15bSPeter Brune     }
805e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
806eed5f15bSPeter Brune     snes->iter = 0;
807eed5f15bSPeter Brune     snes->norm = fnorm;
808e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
809eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
810eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
811eed5f15bSPeter Brune 
812eed5f15bSPeter Brune     /* test convergence */
813eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
814eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
815eed5f15bSPeter Brune   } else {
816e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
817eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
818eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
819eed5f15bSPeter Brune   }
820eed5f15bSPeter Brune 
821eed5f15bSPeter Brune   /* Call general purpose update function */
822eed5f15bSPeter Brune   if (snes->ops->update) {
823eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
824eed5f15bSPeter Brune   }
825eed5f15bSPeter Brune 
826eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
827b22975d2SPatrick Farrell     /* Copy the state before modification by application of the composite solver;
828b22975d2SPatrick Farrell        we will subtract the new state after application */
829b22975d2SPatrick Farrell     ierr = VecCopy(X, snes->work[0]);CHKERRQ(ierr);
830b22975d2SPatrick Farrell 
831eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
832eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
833eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
834eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
83590a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
83690a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
837eed5f15bSPeter Brune     } else {
83890a8ba9bSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
839eed5f15bSPeter Brune     }
840d5a53f06SPeter Brune     if (snes->reason < 0) break;
841d5a53f06SPeter Brune 
842b22975d2SPatrick Farrell     /* Compute the solution update for convergence testing */
843b22975d2SPatrick Farrell     ierr = VecAXPY(snes->work[0], -1.0, X);CHKERRQ(ierr);
844b22975d2SPatrick Farrell     ierr = VecScale(snes->work[0], -1.0);CHKERRQ(ierr);
845b22975d2SPatrick Farrell 
846eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
847eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
848eed5f15bSPeter Brune       if (snes->domainerror) {
849eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
850eed5f15bSPeter Brune         break;
851eed5f15bSPeter Brune       }
852b22975d2SPatrick Farrell 
853a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
854a6da83ecSPatrick Farrell         ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
855a6da83ecSPatrick Farrell         ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
856a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr);
857a6da83ecSPatrick Farrell         ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
858a6da83ecSPatrick Farrell         ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
859a6da83ecSPatrick Farrell       } else {
860b22975d2SPatrick Farrell         ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr);
861b22975d2SPatrick Farrell         ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
862b22975d2SPatrick Farrell         ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
863b22975d2SPatrick Farrell 
864b22975d2SPatrick Farrell         ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr);
865b22975d2SPatrick Farrell         ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
866b22975d2SPatrick Farrell         ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
867a6da83ecSPatrick Farrell       }
868b22975d2SPatrick Farrell 
869eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
870eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
871eed5f15bSPeter Brune         break;
872eed5f15bSPeter Brune       }
873b22975d2SPatrick Farrell     } else if (normtype == SNES_NORM_ALWAYS) {
874b22975d2SPatrick Farrell       ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
875b22975d2SPatrick Farrell       ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
876b22975d2SPatrick Farrell       ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
877b22975d2SPatrick Farrell       ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
878eed5f15bSPeter Brune     }
879eed5f15bSPeter Brune     /* Monitor convergence */
880e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
881eed5f15bSPeter Brune     snes->iter = i+1;
882eed5f15bSPeter Brune     snes->norm = fnorm;
883e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
884eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
885eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
886eed5f15bSPeter Brune     /* Test for convergence */
887b22975d2SPatrick Farrell     if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
888eed5f15bSPeter Brune     if (snes->reason) break;
889eed5f15bSPeter Brune     /* Call general purpose update function */
890eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
891eed5f15bSPeter Brune   }
89272edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS) {
893eed5f15bSPeter Brune     if (i == snes->max_its) {
894eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
895eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
896eed5f15bSPeter Brune     }
897eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
898eed5f15bSPeter Brune   PetscFunctionReturn(0);
899eed5f15bSPeter Brune }
900eed5f15bSPeter Brune 
901eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
902eed5f15bSPeter Brune 
903eed5f15bSPeter Brune /*MC
904eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
905eed5f15bSPeter Brune 
906eed5f15bSPeter Brune    Options Database Keys:
907eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
908eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
909eed5f15bSPeter Brune 
910eed5f15bSPeter Brune    Level: intermediate
911eed5f15bSPeter Brune 
912eed5f15bSPeter Brune    Concepts: composing solvers
913eed5f15bSPeter Brune 
914eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
915eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
916eed5f15bSPeter Brune            SNESCompositeGetSNES()
917eed5f15bSPeter Brune 
918eed5f15bSPeter Brune M*/
919eed5f15bSPeter Brune 
920eed5f15bSPeter Brune #undef __FUNCT__
921eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite"
922eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
923eed5f15bSPeter Brune {
924eed5f15bSPeter Brune   PetscErrorCode ierr;
925eed5f15bSPeter Brune   SNES_Composite   *jac;
926eed5f15bSPeter Brune 
927eed5f15bSPeter Brune   PetscFunctionBegin;
928b00a9115SJed Brown   ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr);
929eed5f15bSPeter Brune 
930eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
931eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
932eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
933eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
934eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
935eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
936eed5f15bSPeter Brune 
937eed5f15bSPeter Brune   snes->data = (void*)jac;
93890a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
93990a8ba9bSPeter Brune   jac->Fes   = NULL;
94090a8ba9bSPeter Brune   jac->Xes   = NULL;
9419c2f3473SPeter Brune   jac->fnorms = NULL;
94290a8ba9bSPeter Brune   jac->nsnes = 0;
943eed5f15bSPeter Brune   jac->head  = 0;
9445e806d2eSPeter Brune   jac->stol  = 0.1;
9455e806d2eSPeter Brune   jac->rtol  = 1.1;
946eed5f15bSPeter Brune 
94790a8ba9bSPeter Brune   jac->h     = NULL;
94890a8ba9bSPeter Brune   jac->s     = NULL;
94990a8ba9bSPeter Brune   jac->beta  = NULL;
95090a8ba9bSPeter Brune   jac->work  = NULL;
95190a8ba9bSPeter Brune   jac->rwork = NULL;
95290a8ba9bSPeter Brune 
9538f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
9548f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
9558f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
9568f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
957eed5f15bSPeter Brune   PetscFunctionReturn(0);
958eed5f15bSPeter Brune }
959eed5f15bSPeter Brune 
960