xref: /petsc/src/snes/impls/composite/snescomposite.c (revision a6da83ec440a2a73c13fe3c1d835dcd2ca29b01d)
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) {
8971dbe336SPeter Brune       ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
90d5a53f06SPeter Brune       if (PetscIsInfOrNanReal(*fnorm)) {
91d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
92d5a53f06SPeter Brune         PetscFunctionReturn(0);
93d5a53f06SPeter Brune       }
94d5a53f06SPeter Brune     }
9572edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
96be95d8f1SBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
97d5a53f06SPeter Brune     if (snes->domainerror) {
98d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
99d5a53f06SPeter Brune       PetscFunctionReturn(0);
100d5a53f06SPeter Brune     }
101d5a53f06SPeter Brune     if (fnorm) {
102*a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
103*a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
104*a6da83ecSPatrick Farrell       } else {
105d5a53f06SPeter Brune         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
106*a6da83ecSPatrick Farrell       }
107*a6da83ecSPatrick Farrell 
108d5a53f06SPeter Brune       if (PetscIsInfOrNanReal(*fnorm)) {
109d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
110d5a53f06SPeter Brune         PetscFunctionReturn(0);
111d5a53f06SPeter Brune       }
112d5a53f06SPeter Brune     }
113eed5f15bSPeter Brune   }
114eed5f15bSPeter Brune   PetscFunctionReturn(0);
115eed5f15bSPeter Brune }
116eed5f15bSPeter Brune 
117eed5f15bSPeter Brune #undef __FUNCT__
118eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive"
119eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
120eed5f15bSPeter Brune {
121eed5f15bSPeter Brune   PetscErrorCode      ierr;
122eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
123eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
124eed5f15bSPeter Brune   Vec                 Y,Xorig;
125d5a53f06SPeter Brune   SNESConvergedReason reason;
126eed5f15bSPeter Brune 
127eed5f15bSPeter Brune   PetscFunctionBegin;
128eed5f15bSPeter Brune   Y = snes->vec_sol_update;
129eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
130eed5f15bSPeter Brune   Xorig = jac->Xorig;
131302440fdSBarry Smith   ierr = VecCopy(X,Xorig);CHKERRQ(ierr);
132eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
13372edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
134eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
135eed5f15bSPeter Brune     while (next->next) {
136eed5f15bSPeter Brune       next = next->next;
137eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
138eed5f15bSPeter Brune     }
139eed5f15bSPeter Brune   }
140eed5f15bSPeter Brune   next = jac->head;
1418f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
142eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
143d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
144d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
145d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
146d5a53f06SPeter Brune     PetscFunctionReturn(0);
147d5a53f06SPeter Brune   }
148eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1498f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1508f626970SPeter Brune   while (next->next) {
1518f626970SPeter Brune     next = next->next;
1528f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1538f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
154d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
155d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
156d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
157d5a53f06SPeter Brune       PetscFunctionReturn(0);
158d5a53f06SPeter Brune     }
1598f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1608f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
161eed5f15bSPeter Brune   }
16272edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
163eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
164*a6da83ecSPatrick Farrell     if (fnorm) {
165*a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
166*a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
167*a6da83ecSPatrick Farrell       } else {
168*a6da83ecSPatrick Farrell         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
169*a6da83ecSPatrick Farrell       }
170*a6da83ecSPatrick Farrell 
171*a6da83ecSPatrick Farrell       if (PetscIsInfOrNanReal(*fnorm)) {
172*a6da83ecSPatrick Farrell         snes->reason = SNES_DIVERGED_FNORM_NAN;
173*a6da83ecSPatrick Farrell         PetscFunctionReturn(0);
174*a6da83ecSPatrick Farrell       }
175*a6da83ecSPatrick Farrell     }
176eed5f15bSPeter Brune   }
177eed5f15bSPeter Brune   PetscFunctionReturn(0);
178eed5f15bSPeter Brune }
17990a8ba9bSPeter Brune 
18090a8ba9bSPeter Brune #undef __FUNCT__
18190a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal"
18290a8ba9bSPeter Brune /* approximately solve the overdetermined system:
18390a8ba9bSPeter Brune 
18490a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
18590a8ba9bSPeter Brune  \alpha_i                      = 1
18690a8ba9bSPeter Brune 
18790a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
18890a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
18990a8ba9bSPeter Brune 
19090a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
19190a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
19290a8ba9bSPeter Brune  */
19390a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
19490a8ba9bSPeter Brune {
19590a8ba9bSPeter Brune   PetscErrorCode      ierr;
19690a8ba9bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
19790a8ba9bSPeter Brune   SNES_CompositeLink  next = jac->head;
19890a8ba9bSPeter Brune   Vec                 *Xes = jac->Xes,*Fes = jac->Fes;
19990a8ba9bSPeter Brune   PetscInt            i,j;
2005e806d2eSPeter Brune   PetscScalar         tot,total,ftf;
2015e806d2eSPeter Brune   PetscReal           min_fnorm;
2025e806d2eSPeter Brune   PetscInt            min_i;
203d5a53f06SPeter Brune   SNESConvergedReason reason;
20490a8ba9bSPeter Brune 
20590a8ba9bSPeter Brune   PetscFunctionBegin;
20690a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
20758bdfa14SPeter Brune 
20858bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
20958bdfa14SPeter Brune     next = jac->head;
21058bdfa14SPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
21158bdfa14SPeter Brune     while (next->next) {
21258bdfa14SPeter Brune       next = next->next;
21358bdfa14SPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
21458bdfa14SPeter Brune     }
21558bdfa14SPeter Brune   }
21658bdfa14SPeter Brune 
21790a8ba9bSPeter Brune   next = jac->head;
21890a8ba9bSPeter Brune   i = 0;
21990a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
22090a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
221d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
222d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
223d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
224d5a53f06SPeter Brune     PetscFunctionReturn(0);
225d5a53f06SPeter Brune   }
22690a8ba9bSPeter Brune   while (next->next) {
22790a8ba9bSPeter Brune     i++;
22890a8ba9bSPeter Brune     next = next->next;
22990a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
23090a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
231d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
232d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
233d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
234d5a53f06SPeter Brune       PetscFunctionReturn(0);
235d5a53f06SPeter Brune     }
23690a8ba9bSPeter Brune   }
23790a8ba9bSPeter Brune 
23890a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
23990a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
24090a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2419c2f3473SPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
24290a8ba9bSPeter Brune     }
2439c2f3473SPeter Brune     ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
24490a8ba9bSPeter Brune   }
2455e806d2eSPeter Brune 
24690a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
24790a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2489c2f3473SPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
2499c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
25090a8ba9bSPeter Brune     }
2519c2f3473SPeter Brune     ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
25290a8ba9bSPeter Brune   }
25390a8ba9bSPeter Brune 
2549c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
25590a8ba9bSPeter Brune 
25690a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
25790a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
2589c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
25990a8ba9bSPeter Brune     }
26090a8ba9bSPeter Brune   }
26190a8ba9bSPeter Brune 
26290a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
2639c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
2649c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
26590a8ba9bSPeter Brune     }
2669c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2679c2f3473SPeter Brune   }
26890a8ba9bSPeter Brune 
26990a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
27090a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
27190a8ba9bSPeter Brune #else
27290a8ba9bSPeter Brune   jac->info  = 0;
27390a8ba9bSPeter Brune   jac->rcond = -1.;
27490a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
27590a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
2769c2f3473SPeter 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));
27790a8ba9bSPeter Brune #else
2789c2f3473SPeter 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));
27990a8ba9bSPeter Brune #endif
28090a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
28190a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
28290a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
28390a8ba9bSPeter Brune #endif
2849c2f3473SPeter Brune   tot = 0.;
2855e806d2eSPeter Brune   total = 0.;
28690a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
28790a8ba9bSPeter Brune     if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
288c783eb9eSBarry Smith     ierr = PetscInfo2(snes,"%D: %g\n",i,(double)PetscRealPart(jac->beta[i]));CHKERRQ(ierr);
2899c2f3473SPeter Brune     tot += jac->beta[i];
2905e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
29190a8ba9bSPeter Brune   }
2929c2f3473SPeter Brune   ierr = VecScale(X,(1. - tot));CHKERRQ(ierr);
29390a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
29490a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
295*a6da83ecSPatrick Farrell 
296*a6da83ecSPatrick Farrell   if (snes->xl && snes->xu) {
297*a6da83ecSPatrick Farrell     ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
298*a6da83ecSPatrick Farrell   } else {
2999c2f3473SPeter Brune     ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
300*a6da83ecSPatrick Farrell   }
30190a8ba9bSPeter Brune 
3025e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
3035e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
3045e806d2eSPeter Brune   min_i     = 0;
3059c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
3065e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
3075e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
3085e806d2eSPeter Brune       min_i     = i;
3099c2f3473SPeter Brune     }
3109c2f3473SPeter Brune   }
3115e806d2eSPeter Brune 
3125e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
3135e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
31458bdfa14SPeter Brune     ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr);
31558bdfa14SPeter Brune     ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr);
3165e806d2eSPeter Brune     *fnorm = min_fnorm;
3175e806d2eSPeter Brune   }
31890a8ba9bSPeter Brune   PetscFunctionReturn(0);
31990a8ba9bSPeter Brune }
32090a8ba9bSPeter Brune 
321eed5f15bSPeter Brune #undef __FUNCT__
322eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite"
323eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
324eed5f15bSPeter Brune {
325eed5f15bSPeter Brune   PetscErrorCode     ierr;
326dd515a93SPeter Brune   DM                 dm;
327eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
328eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
32990a8ba9bSPeter Brune   PetscInt           n=0,i;
33090a8ba9bSPeter Brune   Vec                F;
331eed5f15bSPeter Brune 
332eed5f15bSPeter Brune   PetscFunctionBegin;
333dd515a93SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
334eed5f15bSPeter Brune   while (next) {
33590a8ba9bSPeter Brune     n++;
336dd515a93SPeter Brune     ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr);
337eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
338eed5f15bSPeter Brune     next = next->next;
339eed5f15bSPeter Brune   }
34090a8ba9bSPeter Brune   jac->nsnes = n;
34190a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
34290a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
34390a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
344854ce69bSBarry Smith     ierr = PetscMalloc1(n,&jac->Fes);CHKERRQ(ierr);
345854ce69bSBarry Smith     ierr = PetscMalloc1(n,&jac->fnorms);CHKERRQ(ierr);
34690a8ba9bSPeter Brune     next = jac->head;
34790a8ba9bSPeter Brune     i = 0;
34890a8ba9bSPeter Brune     while (next) {
34990a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
35090a8ba9bSPeter Brune       jac->Fes[i] = F;
35190a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
35290a8ba9bSPeter Brune       next = next->next;
35390a8ba9bSPeter Brune       i++;
35490a8ba9bSPeter Brune     }
35590a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
35690a8ba9bSPeter Brune     jac->nrhs  = 1;
3579c2f3473SPeter Brune     jac->lda   = jac->nsnes;
35890a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
35990a8ba9bSPeter Brune     jac->n     = jac->nsnes;
36090a8ba9bSPeter Brune 
361785e854fSJed Brown     ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr);
362785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr);
363785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr);
364785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr);
3659c2f3473SPeter Brune     jac->lwork = 12*jac->n;
36690a8ba9bSPeter Brune #if PETSC_USE_COMPLEX
367854ce69bSBarry Smith     ierr = PetscMalloc1(jac->lwork,&jac->rwork);CHKERRQ(ierr);
36890a8ba9bSPeter Brune #endif
369854ce69bSBarry Smith     ierr = PetscMalloc1(jac->lwork,&jac->work);CHKERRQ(ierr);
37090a8ba9bSPeter Brune   }
37190a8ba9bSPeter Brune 
372eed5f15bSPeter Brune   PetscFunctionReturn(0);
373eed5f15bSPeter Brune }
374eed5f15bSPeter Brune 
375eed5f15bSPeter Brune #undef __FUNCT__
376eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite"
377eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
378eed5f15bSPeter Brune {
379eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
380eed5f15bSPeter Brune   PetscErrorCode   ierr;
381eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
382eed5f15bSPeter Brune 
383eed5f15bSPeter Brune   PetscFunctionBegin;
384eed5f15bSPeter Brune   while (next) {
385eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
386eed5f15bSPeter Brune     next = next->next;
387eed5f15bSPeter Brune   }
388eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
38990a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
39090a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
3919c2f3473SPeter Brune   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
39290a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
39390a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
3949c2f3473SPeter Brune   ierr = PetscFree(jac->g);CHKERRQ(ierr);
39590a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
39690a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
39790a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
398eed5f15bSPeter Brune   PetscFunctionReturn(0);
399eed5f15bSPeter Brune }
400eed5f15bSPeter Brune 
401eed5f15bSPeter Brune #undef __FUNCT__
402eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite"
403eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
404eed5f15bSPeter Brune {
405eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
406eed5f15bSPeter Brune   PetscErrorCode   ierr;
407eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
408eed5f15bSPeter Brune 
409eed5f15bSPeter Brune   PetscFunctionBegin;
410eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
411eed5f15bSPeter Brune   while (next) {
412eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
413eed5f15bSPeter Brune     next_tmp = next;
414eed5f15bSPeter Brune     next     = next->next;
415eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
416eed5f15bSPeter Brune   }
417eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
418eed5f15bSPeter Brune   PetscFunctionReturn(0);
419eed5f15bSPeter Brune }
420eed5f15bSPeter Brune 
421eed5f15bSPeter Brune #undef __FUNCT__
422eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite"
4238c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(PetscOptions *PetscOptionsObject,SNES snes)
424eed5f15bSPeter Brune {
425eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
426eed5f15bSPeter Brune   PetscErrorCode     ierr;
427eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
428eed5f15bSPeter Brune   SNES_CompositeLink next;
429eed5f15bSPeter Brune   char               *sneses[8];
4308f626970SPeter Brune   PetscReal          dmps[8];
431eed5f15bSPeter Brune   PetscBool          flg;
432eed5f15bSPeter Brune 
433eed5f15bSPeter Brune   PetscFunctionBegin;
434e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");CHKERRQ(ierr);
435eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
436eed5f15bSPeter Brune   if (flg) {
437eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
438eed5f15bSPeter Brune   }
439eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
440eed5f15bSPeter Brune   if (flg) {
441eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
442eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
443eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
444eed5f15bSPeter Brune     }
445eed5f15bSPeter Brune   }
4468f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
4478f626970SPeter Brune   if (flg) {
4488f626970SPeter Brune     for (i=0; i<nmax; i++) {
4498f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
4508f626970SPeter Brune     }
4518f626970SPeter Brune   }
4525e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr);
4535e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr);
454eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
455eed5f15bSPeter Brune 
456eed5f15bSPeter Brune   next = jac->head;
457eed5f15bSPeter Brune   while (next) {
458eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
459eed5f15bSPeter Brune     next = next->next;
460eed5f15bSPeter Brune   }
461eed5f15bSPeter Brune   PetscFunctionReturn(0);
462eed5f15bSPeter Brune }
463eed5f15bSPeter Brune 
464eed5f15bSPeter Brune #undef __FUNCT__
465eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite"
466eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
467eed5f15bSPeter Brune {
468eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
469eed5f15bSPeter Brune   PetscErrorCode   ierr;
470eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
471eed5f15bSPeter Brune   PetscBool        iascii;
472eed5f15bSPeter Brune 
473eed5f15bSPeter Brune   PetscFunctionBegin;
474eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
475eed5f15bSPeter Brune   if (iascii) {
476eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
477eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
478eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
479eed5f15bSPeter Brune   }
480eed5f15bSPeter Brune   if (iascii) {
481eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
482eed5f15bSPeter Brune   }
483eed5f15bSPeter Brune   while (next) {
484eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
485eed5f15bSPeter Brune     next = next->next;
486eed5f15bSPeter Brune   }
487eed5f15bSPeter Brune   if (iascii) {
488eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
489eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
490eed5f15bSPeter Brune   }
491eed5f15bSPeter Brune   PetscFunctionReturn(0);
492eed5f15bSPeter Brune }
493eed5f15bSPeter Brune 
494eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
495eed5f15bSPeter Brune 
496eed5f15bSPeter Brune #undef __FUNCT__
497eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite"
498eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
499eed5f15bSPeter Brune {
500eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
501eed5f15bSPeter Brune 
502eed5f15bSPeter Brune   PetscFunctionBegin;
503eed5f15bSPeter Brune   jac->type = type;
504eed5f15bSPeter Brune   PetscFunctionReturn(0);
505eed5f15bSPeter Brune }
506eed5f15bSPeter Brune 
507eed5f15bSPeter Brune #undef __FUNCT__
508eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite"
509eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
510eed5f15bSPeter Brune {
511eed5f15bSPeter Brune   SNES_Composite     *jac;
512eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
513eed5f15bSPeter Brune   PetscErrorCode   ierr;
514eed5f15bSPeter Brune   PetscInt         cnt = 0;
515eed5f15bSPeter Brune   const char       *prefix;
516eed5f15bSPeter Brune   char             newprefix[8];
517eed5f15bSPeter Brune   DM               dm;
518eed5f15bSPeter Brune 
519eed5f15bSPeter Brune   PetscFunctionBegin;
520b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ilink);CHKERRQ(ierr);
521eed5f15bSPeter Brune   ilink->next = 0;
522eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
523eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
524eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
525eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
526cf5b3eb5SPeter Brune   ierr        = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr);
527eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
528eed5f15bSPeter Brune   next = jac->head;
529eed5f15bSPeter Brune   if (!next) {
530eed5f15bSPeter Brune     jac->head       = ilink;
531eed5f15bSPeter Brune     ilink->previous = NULL;
532eed5f15bSPeter Brune   } else {
533eed5f15bSPeter Brune     cnt++;
534eed5f15bSPeter Brune     while (next->next) {
535eed5f15bSPeter Brune       next = next->next;
536eed5f15bSPeter Brune       cnt++;
537eed5f15bSPeter Brune     }
538eed5f15bSPeter Brune     next->next      = ilink;
539eed5f15bSPeter Brune     ilink->previous = next;
540eed5f15bSPeter Brune   }
541eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
542eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
543eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
544eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
5458f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
546eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
54772edecb9SPeter Brune   ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
5488f626970SPeter Brune   ilink->dmp = 1.0;
54990a8ba9bSPeter Brune   jac->nsnes++;
550eed5f15bSPeter Brune   PetscFunctionReturn(0);
551eed5f15bSPeter Brune }
552eed5f15bSPeter Brune 
553eed5f15bSPeter Brune #undef __FUNCT__
554eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite"
555eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
556eed5f15bSPeter Brune {
557eed5f15bSPeter Brune   SNES_Composite     *jac;
558eed5f15bSPeter Brune   SNES_CompositeLink next;
559eed5f15bSPeter Brune   PetscInt         i;
560eed5f15bSPeter Brune 
561eed5f15bSPeter Brune   PetscFunctionBegin;
562eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
563eed5f15bSPeter Brune   next = jac->head;
564eed5f15bSPeter Brune   for (i=0; i<n; i++) {
565eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
566eed5f15bSPeter Brune     next = next->next;
567eed5f15bSPeter Brune   }
568eed5f15bSPeter Brune   *subsnes = next->snes;
569eed5f15bSPeter Brune   PetscFunctionReturn(0);
570eed5f15bSPeter Brune }
571eed5f15bSPeter Brune 
572eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
573eed5f15bSPeter Brune #undef __FUNCT__
574eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType"
575e31ab4f9SPeter Brune /*@C
576eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
577eed5f15bSPeter Brune 
578eed5f15bSPeter Brune    Logically Collective on SNES
579eed5f15bSPeter Brune 
580eed5f15bSPeter Brune    Input Parameter:
581eed5f15bSPeter Brune +  snes - the preconditioner context
582eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
583eed5f15bSPeter Brune 
584eed5f15bSPeter Brune    Options Database Key:
585eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
586eed5f15bSPeter Brune 
587eed5f15bSPeter Brune    Level: Developer
588eed5f15bSPeter Brune 
589eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
590eed5f15bSPeter Brune @*/
591eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
592eed5f15bSPeter Brune {
593eed5f15bSPeter Brune   PetscErrorCode ierr;
594eed5f15bSPeter Brune 
595eed5f15bSPeter Brune   PetscFunctionBegin;
596eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
597eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
598eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
599eed5f15bSPeter Brune   PetscFunctionReturn(0);
600eed5f15bSPeter Brune }
601eed5f15bSPeter Brune 
602eed5f15bSPeter Brune #undef __FUNCT__
603eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES"
604eed5f15bSPeter Brune /*@C
605eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
606eed5f15bSPeter Brune 
607eed5f15bSPeter Brune    Collective on SNES
608eed5f15bSPeter Brune 
609eed5f15bSPeter Brune    Input Parameters:
610eed5f15bSPeter Brune +  snes - the preconditioner context
611eed5f15bSPeter Brune -  type - the type of the new preconditioner
612eed5f15bSPeter Brune 
613eed5f15bSPeter Brune    Level: Developer
614eed5f15bSPeter Brune 
615eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add
616eed5f15bSPeter Brune @*/
617eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
618eed5f15bSPeter Brune {
619eed5f15bSPeter Brune   PetscErrorCode ierr;
620eed5f15bSPeter Brune 
621eed5f15bSPeter Brune   PetscFunctionBegin;
622eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
623eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
624eed5f15bSPeter Brune   PetscFunctionReturn(0);
625eed5f15bSPeter Brune }
626eed5f15bSPeter Brune #undef __FUNCT__
627eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES"
628eed5f15bSPeter Brune /*@
629eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
630eed5f15bSPeter Brune 
631eed5f15bSPeter Brune    Not Collective
632eed5f15bSPeter Brune 
633eed5f15bSPeter Brune    Input Parameter:
634eed5f15bSPeter Brune +  snes - the preconditioner context
635eed5f15bSPeter Brune -  n - the number of the snes requested
636eed5f15bSPeter Brune 
637eed5f15bSPeter Brune    Output Parameters:
638eed5f15bSPeter Brune .  subsnes - the SNES requested
639eed5f15bSPeter Brune 
640eed5f15bSPeter Brune    Level: Developer
641eed5f15bSPeter Brune 
642eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
643eed5f15bSPeter Brune 
644eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
645eed5f15bSPeter Brune @*/
646eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
647eed5f15bSPeter Brune {
648eed5f15bSPeter Brune   PetscErrorCode ierr;
649eed5f15bSPeter Brune 
650eed5f15bSPeter Brune   PetscFunctionBegin;
651eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
652eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
653eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
654eed5f15bSPeter Brune   PetscFunctionReturn(0);
655eed5f15bSPeter Brune }
656eed5f15bSPeter Brune 
657eed5f15bSPeter Brune #undef __FUNCT__
6586b2b7f7bSPatrick Farrell #define __FUNCT__ "SNESCompositeGetNumber"
6596b2b7f7bSPatrick Farrell /*@
6606b2b7f7bSPatrick Farrell    SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES.
6616b2b7f7bSPatrick Farrell 
6626b2b7f7bSPatrick Farrell    Logically Collective on SNES
6636b2b7f7bSPatrick Farrell 
6646b2b7f7bSPatrick Farrell    Input Parameter:
6656b2b7f7bSPatrick Farrell    snes - the preconditioner context
6666b2b7f7bSPatrick Farrell 
6676b2b7f7bSPatrick Farrell    Output Parameter:
6686b2b7f7bSPatrick Farrell    n - the number of subsolvers
6696b2b7f7bSPatrick Farrell 
6706b2b7f7bSPatrick Farrell    Level: Developer
6716b2b7f7bSPatrick Farrell 
6726b2b7f7bSPatrick Farrell .keywords: SNES, composite preconditioner
6736b2b7f7bSPatrick Farrell @*/
6746b2b7f7bSPatrick Farrell PetscErrorCode  SNESCompositeGetNumber(SNES snes,PetscInt *n)
6756b2b7f7bSPatrick Farrell {
6766b2b7f7bSPatrick Farrell   SNES_Composite     *jac;
6776b2b7f7bSPatrick Farrell   SNES_CompositeLink next;
6786b2b7f7bSPatrick Farrell 
6796b2b7f7bSPatrick Farrell   PetscFunctionBegin;
6806b2b7f7bSPatrick Farrell   jac  = (SNES_Composite*)snes->data;
6816b2b7f7bSPatrick Farrell   next = jac->head;
6826b2b7f7bSPatrick Farrell 
6836b2b7f7bSPatrick Farrell   *n = 0;
6846b2b7f7bSPatrick Farrell   while (next) {
6856b2b7f7bSPatrick Farrell     *n = *n + 1;
6866b2b7f7bSPatrick Farrell     next = next->next;
6876b2b7f7bSPatrick Farrell   }
6886b2b7f7bSPatrick Farrell   PetscFunctionReturn(0);
6896b2b7f7bSPatrick Farrell }
6906b2b7f7bSPatrick Farrell 
6916b2b7f7bSPatrick Farrell #undef __FUNCT__
6928f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite"
6938f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
6948f626970SPeter Brune {
6958f626970SPeter Brune   SNES_Composite     *jac;
6968f626970SPeter Brune   SNES_CompositeLink next;
6978f626970SPeter Brune   PetscInt         i;
6988f626970SPeter Brune 
6998f626970SPeter Brune   PetscFunctionBegin;
7008f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
7018f626970SPeter Brune   next = jac->head;
7028f626970SPeter Brune   for (i=0; i<n; i++) {
7038f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
7048f626970SPeter Brune     next = next->next;
7058f626970SPeter Brune   }
7068f626970SPeter Brune   next->dmp = dmp;
7078f626970SPeter Brune   PetscFunctionReturn(0);
7088f626970SPeter Brune }
7098f626970SPeter Brune 
7108f626970SPeter Brune #undef __FUNCT__
7118f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping"
7128f626970SPeter Brune /*@
7138f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
7148f626970SPeter Brune 
7158f626970SPeter Brune    Not Collective
7168f626970SPeter Brune 
7178f626970SPeter Brune    Input Parameter:
7188f626970SPeter Brune +  snes - the preconditioner context
7198f626970SPeter Brune .  n - the number of the snes requested
7208f626970SPeter Brune -  dmp - the damping
7218f626970SPeter Brune 
7228f626970SPeter Brune    Level: Developer
7238f626970SPeter Brune 
7248f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
7258f626970SPeter Brune 
7268f626970SPeter Brune .seealso: SNESCompositeAddSNES()
7278f626970SPeter Brune @*/
7288f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
7298f626970SPeter Brune {
7308f626970SPeter Brune   PetscErrorCode ierr;
7318f626970SPeter Brune 
7328f626970SPeter Brune   PetscFunctionBegin;
7338f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7348f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
7358f626970SPeter Brune   PetscFunctionReturn(0);
7368f626970SPeter Brune }
7378f626970SPeter Brune 
7388f626970SPeter Brune #undef __FUNCT__
739eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite"
740eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes)
741eed5f15bSPeter Brune {
742eed5f15bSPeter Brune   Vec            F;
743eed5f15bSPeter Brune   Vec            X;
744eed5f15bSPeter Brune   Vec            B;
745eed5f15bSPeter Brune   PetscInt       i;
746b22975d2SPatrick Farrell   PetscReal      fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
747eed5f15bSPeter Brune   PetscErrorCode ierr;
74872edecb9SPeter Brune   SNESNormSchedule normtype;
749eed5f15bSPeter Brune   SNES_Composite *comp = (SNES_Composite*)snes->data;
750eed5f15bSPeter Brune 
751eed5f15bSPeter Brune   PetscFunctionBegin;
752eed5f15bSPeter Brune   X = snes->vec_sol;
753eed5f15bSPeter Brune   F = snes->vec_func;
754eed5f15bSPeter Brune   B = snes->vec_rhs;
755eed5f15bSPeter Brune 
756e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
757eed5f15bSPeter Brune   snes->iter   = 0;
758eed5f15bSPeter Brune   snes->norm   = 0.;
759e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
760b22975d2SPatrick Farrell   ierr         = SNESSetWorkVecs(snes, 1);CHKERRQ(ierr);
761eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
76272edecb9SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr);
76372edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
764eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
765eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
766eed5f15bSPeter Brune       if (snes->domainerror) {
767eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
768eed5f15bSPeter Brune         PetscFunctionReturn(0);
769eed5f15bSPeter Brune       }
770eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
771eed5f15bSPeter Brune 
772*a6da83ecSPatrick Farrell     if (snes->xl && snes->xu) {
773*a6da83ecSPatrick Farrell       ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr);
774*a6da83ecSPatrick Farrell     } else {
775eed5f15bSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
776*a6da83ecSPatrick Farrell     }
777*a6da83ecSPatrick Farrell 
778eed5f15bSPeter Brune     if (PetscIsInfOrNanReal(fnorm)) {
779eed5f15bSPeter Brune       snes->reason = SNES_DIVERGED_FNORM_NAN;
780eed5f15bSPeter Brune       PetscFunctionReturn(0);
781eed5f15bSPeter Brune     }
782e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
783eed5f15bSPeter Brune     snes->iter = 0;
784eed5f15bSPeter Brune     snes->norm = fnorm;
785e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
786eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
787eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
788eed5f15bSPeter Brune 
789eed5f15bSPeter Brune     /* test convergence */
790eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
791eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
792eed5f15bSPeter Brune   } else {
793e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
794eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
795eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
796eed5f15bSPeter Brune   }
797eed5f15bSPeter Brune 
798eed5f15bSPeter Brune   /* Call general purpose update function */
799eed5f15bSPeter Brune   if (snes->ops->update) {
800eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
801eed5f15bSPeter Brune   }
802eed5f15bSPeter Brune 
803eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
804b22975d2SPatrick Farrell     /* Copy the state before modification by application of the composite solver;
805b22975d2SPatrick Farrell        we will subtract the new state after application */
806b22975d2SPatrick Farrell     ierr = VecCopy(X, snes->work[0]);CHKERRQ(ierr);
807b22975d2SPatrick Farrell 
808eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
809eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
810eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
811eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
81290a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
81390a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
814eed5f15bSPeter Brune     } else {
81590a8ba9bSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
816eed5f15bSPeter Brune     }
817d5a53f06SPeter Brune     if (snes->reason < 0) break;
818d5a53f06SPeter Brune 
819b22975d2SPatrick Farrell     /* Compute the solution update for convergence testing */
820b22975d2SPatrick Farrell     ierr = VecAXPY(snes->work[0], -1.0, X);CHKERRQ(ierr);
821b22975d2SPatrick Farrell     ierr = VecScale(snes->work[0], -1.0);CHKERRQ(ierr);
822b22975d2SPatrick Farrell 
823eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
824eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
825eed5f15bSPeter Brune       if (snes->domainerror) {
826eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
827eed5f15bSPeter Brune         break;
828eed5f15bSPeter Brune       }
829b22975d2SPatrick Farrell 
830*a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
831*a6da83ecSPatrick Farrell         ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
832*a6da83ecSPatrick Farrell         ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
833*a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr);
834*a6da83ecSPatrick Farrell         ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
835*a6da83ecSPatrick Farrell         ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
836*a6da83ecSPatrick Farrell       } else {
837b22975d2SPatrick Farrell         ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr);
838b22975d2SPatrick Farrell         ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
839b22975d2SPatrick Farrell         ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
840b22975d2SPatrick Farrell 
841b22975d2SPatrick Farrell         ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr);
842b22975d2SPatrick Farrell         ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
843b22975d2SPatrick Farrell         ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
844*a6da83ecSPatrick Farrell       }
845b22975d2SPatrick Farrell 
846eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
847eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
848eed5f15bSPeter Brune         break;
849eed5f15bSPeter Brune       }
850b22975d2SPatrick Farrell     } else if (normtype == SNES_NORM_ALWAYS) {
851b22975d2SPatrick Farrell       ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
852b22975d2SPatrick Farrell       ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
853b22975d2SPatrick Farrell       ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
854b22975d2SPatrick Farrell       ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
855eed5f15bSPeter Brune     }
856eed5f15bSPeter Brune     /* Monitor convergence */
857e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
858eed5f15bSPeter Brune     snes->iter = i+1;
859eed5f15bSPeter Brune     snes->norm = fnorm;
860e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
861eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
862eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
863eed5f15bSPeter Brune     /* Test for convergence */
864b22975d2SPatrick Farrell     if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
865eed5f15bSPeter Brune     if (snes->reason) break;
866eed5f15bSPeter Brune     /* Call general purpose update function */
867eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
868eed5f15bSPeter Brune   }
86972edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS) {
870eed5f15bSPeter Brune     if (i == snes->max_its) {
871eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
872eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
873eed5f15bSPeter Brune     }
874eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
875eed5f15bSPeter Brune   PetscFunctionReturn(0);
876eed5f15bSPeter Brune }
877eed5f15bSPeter Brune 
878eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
879eed5f15bSPeter Brune 
880eed5f15bSPeter Brune /*MC
881eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
882eed5f15bSPeter Brune 
883eed5f15bSPeter Brune    Options Database Keys:
884eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
885eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
886eed5f15bSPeter Brune 
887eed5f15bSPeter Brune    Level: intermediate
888eed5f15bSPeter Brune 
889eed5f15bSPeter Brune    Concepts: composing solvers
890eed5f15bSPeter Brune 
891eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
892eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
893eed5f15bSPeter Brune            SNESCompositeGetSNES()
894eed5f15bSPeter Brune 
895eed5f15bSPeter Brune M*/
896eed5f15bSPeter Brune 
897eed5f15bSPeter Brune #undef __FUNCT__
898eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite"
899eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
900eed5f15bSPeter Brune {
901eed5f15bSPeter Brune   PetscErrorCode ierr;
902eed5f15bSPeter Brune   SNES_Composite   *jac;
903eed5f15bSPeter Brune 
904eed5f15bSPeter Brune   PetscFunctionBegin;
905b00a9115SJed Brown   ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr);
906eed5f15bSPeter Brune 
907eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
908eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
909eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
910eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
911eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
912eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
913eed5f15bSPeter Brune 
914eed5f15bSPeter Brune   snes->data = (void*)jac;
91590a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
91690a8ba9bSPeter Brune   jac->Fes   = NULL;
91790a8ba9bSPeter Brune   jac->Xes   = NULL;
9189c2f3473SPeter Brune   jac->fnorms = NULL;
91990a8ba9bSPeter Brune   jac->nsnes = 0;
920eed5f15bSPeter Brune   jac->head  = 0;
9215e806d2eSPeter Brune   jac->stol  = 0.1;
9225e806d2eSPeter Brune   jac->rtol  = 1.1;
923eed5f15bSPeter Brune 
92490a8ba9bSPeter Brune   jac->h     = NULL;
92590a8ba9bSPeter Brune   jac->s     = NULL;
92690a8ba9bSPeter Brune   jac->beta  = NULL;
92790a8ba9bSPeter Brune   jac->work  = NULL;
92890a8ba9bSPeter Brune   jac->rwork = NULL;
92990a8ba9bSPeter Brune 
9308f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
9318f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
9328f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
9338f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
934eed5f15bSPeter Brune   PetscFunctionReturn(0);
935eed5f15bSPeter Brune }
936eed5f15bSPeter Brune 
937