xref: /petsc/src/snes/impls/composite/snescomposite.c (revision c1e67a491dcd3fce18efcb05d6f137bc2aa8e6a9)
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;
240b762f1fSPatrick Farrell   PetscInt           innerFailures; /* the number of inner failures we've seen */
2590a8ba9bSPeter Brune 
2690a8ba9bSPeter Brune   /* context for ADDITIVEOPTIMAL */
2790a8ba9bSPeter Brune   Vec                *Xes,*Fes;      /* solution and residual vectors for the subsolvers */
289c2f3473SPeter Brune   PetscReal          *fnorms;        /* norms of the solutions */
2990a8ba9bSPeter Brune   PetscScalar        *h;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
309c2f3473SPeter Brune   PetscScalar        *g;             /* the dotproducts of the previous function with the candidate functions */
319c2f3473SPeter Brune   PetscBLASInt       n;              /* matrix dimension -- nsnes */
3290a8ba9bSPeter Brune   PetscBLASInt       nrhs;           /* the number of right hand sides */
3390a8ba9bSPeter Brune   PetscBLASInt       lda;            /* the padded matrix dimension */
3490a8ba9bSPeter Brune   PetscBLASInt       ldb;            /* the padded vector dimension */
3590a8ba9bSPeter Brune   PetscReal          *s;             /* the singular values */
369c2f3473SPeter Brune   PetscScalar        *beta;          /* the RHS and combination */
3790a8ba9bSPeter Brune   PetscReal          rcond;          /* the exit condition */
3890a8ba9bSPeter Brune   PetscBLASInt       rank;           /* the effective rank */
3990a8ba9bSPeter Brune   PetscScalar        *work;          /* the work vector */
4090a8ba9bSPeter Brune   PetscReal          *rwork;         /* the real work vector used for complex */
4190a8ba9bSPeter Brune   PetscBLASInt       lwork;          /* the size of the work vector */
4290a8ba9bSPeter Brune   PetscBLASInt       info;           /* the output condition */
4390a8ba9bSPeter Brune 
445e806d2eSPeter Brune   PetscReal          rtol;           /* restart tolerance for accepting the combination */
455e806d2eSPeter Brune   PetscReal          stol;           /* restart tolerance for the combination */
46eed5f15bSPeter Brune } SNES_Composite;
47eed5f15bSPeter Brune 
48eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
49eed5f15bSPeter Brune {
50eed5f15bSPeter Brune   PetscErrorCode      ierr;
51eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
52eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
53eed5f15bSPeter Brune   Vec                 FSub;
54d5a53f06SPeter Brune   SNESConvergedReason reason;
55eed5f15bSPeter Brune 
56eed5f15bSPeter Brune   PetscFunctionBegin;
57eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
5872edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
59eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
60eed5f15bSPeter Brune   }
61eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
62d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
63d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
640b762f1fSPatrick Farrell     jac->innerFailures++;
650b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
66d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
67d5a53f06SPeter Brune       PetscFunctionReturn(0);
68d5a53f06SPeter Brune     }
690b762f1fSPatrick Farrell   }
70eed5f15bSPeter Brune 
71eed5f15bSPeter Brune   while (next->next) {
72eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
73efd4aadfSBarry Smith     if (next->snes->npcside== PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
74eed5f15bSPeter Brune       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
75eed5f15bSPeter Brune       next = next->next;
76eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
77eed5f15bSPeter Brune     } else {
78eed5f15bSPeter Brune       next = next->next;
79eed5f15bSPeter Brune     }
80eed5f15bSPeter Brune     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
81d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
82d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
830b762f1fSPatrick Farrell       jac->innerFailures++;
840b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
85d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
86d5a53f06SPeter Brune         PetscFunctionReturn(0);
87d5a53f06SPeter Brune       }
88eed5f15bSPeter Brune     }
890b762f1fSPatrick Farrell   }
90efd4aadfSBarry Smith   if (next->snes->npcside== PC_RIGHT) {
91eed5f15bSPeter Brune     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
92eed5f15bSPeter Brune     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
93d5a53f06SPeter Brune     if (fnorm) {
9463cdc2bdSPatrick Farrell       if (snes->xl && snes->xu) {
9563cdc2bdSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
9663cdc2bdSPatrick Farrell       } else {
9771dbe336SPeter Brune         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
9863cdc2bdSPatrick Farrell       }
99422a814eSBarry Smith       SNESCheckFunctionNorm(snes,*fnorm);
100d5a53f06SPeter Brune     }
10172edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
102be95d8f1SBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
103d5a53f06SPeter Brune     if (fnorm) {
104a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
105a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
106a6da83ecSPatrick Farrell       } else {
107d5a53f06SPeter Brune         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
108a6da83ecSPatrick Farrell       }
109422a814eSBarry Smith       SNESCheckFunctionNorm(snes,*fnorm);
110d5a53f06SPeter Brune     }
111eed5f15bSPeter Brune   }
112eed5f15bSPeter Brune   PetscFunctionReturn(0);
113eed5f15bSPeter Brune }
114eed5f15bSPeter Brune 
115eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
116eed5f15bSPeter Brune {
117eed5f15bSPeter Brune   PetscErrorCode      ierr;
118eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
119eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
120eed5f15bSPeter Brune   Vec                 Y,Xorig;
121d5a53f06SPeter Brune   SNESConvergedReason reason;
122eed5f15bSPeter Brune 
123eed5f15bSPeter Brune   PetscFunctionBegin;
124eed5f15bSPeter Brune   Y = snes->vec_sol_update;
125eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
126eed5f15bSPeter Brune   Xorig = jac->Xorig;
127302440fdSBarry Smith   ierr = VecCopy(X,Xorig);CHKERRQ(ierr);
128eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
12972edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
130eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
131eed5f15bSPeter Brune     while (next->next) {
132eed5f15bSPeter Brune       next = next->next;
133eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
134eed5f15bSPeter Brune     }
135eed5f15bSPeter Brune   }
136eed5f15bSPeter Brune   next = jac->head;
1378f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
138eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
139d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
140d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1410b762f1fSPatrick Farrell     jac->innerFailures++;
1420b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
143d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
144d5a53f06SPeter Brune       PetscFunctionReturn(0);
145d5a53f06SPeter Brune     }
1460b762f1fSPatrick Farrell   }
147eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1488f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1498f626970SPeter Brune   while (next->next) {
1508f626970SPeter Brune     next = next->next;
1518f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1528f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
153d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
154d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1550b762f1fSPatrick Farrell       jac->innerFailures++;
1560b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
157d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
158d5a53f06SPeter Brune         PetscFunctionReturn(0);
159d5a53f06SPeter Brune       }
1600b762f1fSPatrick Farrell     }
1618f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1628f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
163eed5f15bSPeter Brune   }
16472edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
165eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
166a6da83ecSPatrick Farrell     if (fnorm) {
167a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
168a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
169a6da83ecSPatrick Farrell       } else {
170a6da83ecSPatrick Farrell         ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
171a6da83ecSPatrick Farrell       }
172422a814eSBarry Smith       SNESCheckFunctionNorm(snes,*fnorm);
173a6da83ecSPatrick Farrell     }
174eed5f15bSPeter Brune   }
175eed5f15bSPeter Brune   PetscFunctionReturn(0);
176eed5f15bSPeter Brune }
17790a8ba9bSPeter Brune 
17890a8ba9bSPeter Brune /* approximately solve the overdetermined system:
17990a8ba9bSPeter Brune 
18090a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
18190a8ba9bSPeter Brune  \alpha_i                      = 1
18290a8ba9bSPeter Brune 
18390a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
18490a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
18590a8ba9bSPeter Brune 
18690a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
18790a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
18890a8ba9bSPeter Brune  */
18990a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
19090a8ba9bSPeter Brune {
19190a8ba9bSPeter Brune   PetscErrorCode      ierr;
19290a8ba9bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
19390a8ba9bSPeter Brune   SNES_CompositeLink  next = jac->head;
19490a8ba9bSPeter Brune   Vec                 *Xes = jac->Xes,*Fes = jac->Fes;
19590a8ba9bSPeter Brune   PetscInt            i,j;
1965e806d2eSPeter Brune   PetscScalar         tot,total,ftf;
1975e806d2eSPeter Brune   PetscReal           min_fnorm;
1985e806d2eSPeter Brune   PetscInt            min_i;
199d5a53f06SPeter Brune   SNESConvergedReason reason;
20090a8ba9bSPeter Brune 
20190a8ba9bSPeter Brune   PetscFunctionBegin;
20290a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
20358bdfa14SPeter Brune 
20458bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
20558bdfa14SPeter Brune     next = jac->head;
20658bdfa14SPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
20758bdfa14SPeter Brune     while (next->next) {
20858bdfa14SPeter Brune       next = next->next;
20958bdfa14SPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
21058bdfa14SPeter Brune     }
21158bdfa14SPeter Brune   }
21258bdfa14SPeter Brune 
21390a8ba9bSPeter Brune   next = jac->head;
21490a8ba9bSPeter Brune   i = 0;
21590a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
21690a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
217d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
218d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2190b762f1fSPatrick Farrell     jac->innerFailures++;
2200b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
221d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
222d5a53f06SPeter Brune       PetscFunctionReturn(0);
223d5a53f06SPeter Brune     }
2240b762f1fSPatrick Farrell   }
22590a8ba9bSPeter Brune   while (next->next) {
22690a8ba9bSPeter Brune     i++;
22790a8ba9bSPeter Brune     next = next->next;
22890a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
22990a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
230d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
231d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2320b762f1fSPatrick Farrell       jac->innerFailures++;
2330b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
234d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
235d5a53f06SPeter Brune         PetscFunctionReturn(0);
236d5a53f06SPeter Brune       }
23790a8ba9bSPeter Brune     }
2380b762f1fSPatrick Farrell   }
23990a8ba9bSPeter Brune 
24090a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
24190a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
24290a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2439c2f3473SPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
24490a8ba9bSPeter Brune     }
2459c2f3473SPeter Brune     ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
24690a8ba9bSPeter Brune   }
2475e806d2eSPeter Brune 
24890a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
24990a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2509c2f3473SPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
2519c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
25290a8ba9bSPeter Brune     }
2539c2f3473SPeter Brune     ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
25490a8ba9bSPeter Brune   }
25590a8ba9bSPeter Brune 
2569c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
25790a8ba9bSPeter Brune 
25890a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
25990a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
2609c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
26190a8ba9bSPeter Brune     }
26290a8ba9bSPeter Brune   }
26390a8ba9bSPeter Brune 
26490a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
2659c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
2669c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
26790a8ba9bSPeter Brune     }
2689c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2699c2f3473SPeter Brune   }
27090a8ba9bSPeter Brune 
27190a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
27290a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
27390a8ba9bSPeter Brune #else
27490a8ba9bSPeter Brune   jac->info  = 0;
27590a8ba9bSPeter Brune   jac->rcond = -1.;
27690a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
27790a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
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->rwork,&jac->info));
27990a8ba9bSPeter Brune #else
2809c2f3473SPeter 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));
28190a8ba9bSPeter Brune #endif
28290a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
28390a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
28490a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
28590a8ba9bSPeter Brune #endif
2869c2f3473SPeter Brune   tot = 0.;
2875e806d2eSPeter Brune   total = 0.;
28890a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
289422a814eSBarry Smith     if (snes->errorifnotconverged && PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
290c783eb9eSBarry Smith     ierr = PetscInfo2(snes,"%D: %g\n",i,(double)PetscRealPart(jac->beta[i]));CHKERRQ(ierr);
2919c2f3473SPeter Brune     tot += jac->beta[i];
2925e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
29390a8ba9bSPeter Brune   }
2949c2f3473SPeter Brune   ierr = VecScale(X,(1. - tot));CHKERRQ(ierr);
29590a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
29690a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
297a6da83ecSPatrick Farrell 
298a6da83ecSPatrick Farrell   if (snes->xl && snes->xu) {
299a6da83ecSPatrick Farrell     ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr);
300a6da83ecSPatrick Farrell   } else {
3019c2f3473SPeter Brune     ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
302a6da83ecSPatrick Farrell   }
30390a8ba9bSPeter Brune 
3045e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
3055e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
3065e806d2eSPeter Brune   min_i     = 0;
3079c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
3085e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
3095e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
3105e806d2eSPeter Brune       min_i     = i;
3119c2f3473SPeter Brune     }
3129c2f3473SPeter Brune   }
3135e806d2eSPeter Brune 
3145e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
3155e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
31658bdfa14SPeter Brune     ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr);
31758bdfa14SPeter Brune     ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr);
3185e806d2eSPeter Brune     *fnorm = min_fnorm;
3195e806d2eSPeter Brune   }
32090a8ba9bSPeter Brune   PetscFunctionReturn(0);
32190a8ba9bSPeter Brune }
32290a8ba9bSPeter Brune 
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);
33463cdc2bdSPatrick Farrell 
33563cdc2bdSPatrick Farrell   if (snes->ops->computevariablebounds) {
33663cdc2bdSPatrick Farrell     /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */
33763cdc2bdSPatrick Farrell     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
33863cdc2bdSPatrick Farrell     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
33963cdc2bdSPatrick Farrell     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
34063cdc2bdSPatrick Farrell   }
34163cdc2bdSPatrick Farrell 
342eed5f15bSPeter Brune   while (next) {
34390a8ba9bSPeter Brune     n++;
344dd515a93SPeter Brune     ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr);
34563cdc2bdSPatrick Farrell     ierr = SNESSetApplicationContext(next->snes, snes->user);CHKERRQ(ierr);
34663cdc2bdSPatrick Farrell     if (snes->xl && snes->xu) {
34763cdc2bdSPatrick Farrell       if (snes->ops->computevariablebounds) {
34863cdc2bdSPatrick Farrell         ierr = SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds);CHKERRQ(ierr);
34963cdc2bdSPatrick Farrell       } else {
35063cdc2bdSPatrick Farrell         ierr = SNESVISetVariableBounds(next->snes,snes->xl,snes->xu);CHKERRQ(ierr);
35163cdc2bdSPatrick Farrell       }
35263cdc2bdSPatrick Farrell     }
35363cdc2bdSPatrick Farrell 
354eed5f15bSPeter Brune     next = next->next;
355eed5f15bSPeter Brune   }
35690a8ba9bSPeter Brune   jac->nsnes = n;
35790a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
35890a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
35990a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
360854ce69bSBarry Smith     ierr = PetscMalloc1(n,&jac->Fes);CHKERRQ(ierr);
361854ce69bSBarry Smith     ierr = PetscMalloc1(n,&jac->fnorms);CHKERRQ(ierr);
36290a8ba9bSPeter Brune     next = jac->head;
36390a8ba9bSPeter Brune     i = 0;
36490a8ba9bSPeter Brune     while (next) {
36590a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
36690a8ba9bSPeter Brune       jac->Fes[i] = F;
36790a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
36890a8ba9bSPeter Brune       next = next->next;
36990a8ba9bSPeter Brune       i++;
37090a8ba9bSPeter Brune     }
37190a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
37290a8ba9bSPeter Brune     jac->nrhs  = 1;
3739c2f3473SPeter Brune     jac->lda   = jac->nsnes;
37490a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
37590a8ba9bSPeter Brune     jac->n     = jac->nsnes;
37690a8ba9bSPeter Brune 
377785e854fSJed Brown     ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr);
378785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr);
379785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr);
380785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr);
3819c2f3473SPeter Brune     jac->lwork = 12*jac->n;
382dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
383854ce69bSBarry Smith     ierr = PetscMalloc1(jac->lwork,&jac->rwork);CHKERRQ(ierr);
38490a8ba9bSPeter Brune #endif
385854ce69bSBarry Smith     ierr = PetscMalloc1(jac->lwork,&jac->work);CHKERRQ(ierr);
38690a8ba9bSPeter Brune   }
38790a8ba9bSPeter Brune 
388eed5f15bSPeter Brune   PetscFunctionReturn(0);
389eed5f15bSPeter Brune }
390eed5f15bSPeter Brune 
391eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
392eed5f15bSPeter Brune {
393eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
394eed5f15bSPeter Brune   PetscErrorCode   ierr;
395eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
396eed5f15bSPeter Brune 
397eed5f15bSPeter Brune   PetscFunctionBegin;
398eed5f15bSPeter Brune   while (next) {
399eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
400eed5f15bSPeter Brune     next = next->next;
401eed5f15bSPeter Brune   }
402eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
40390a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
40490a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
4059c2f3473SPeter Brune   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
40690a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
40790a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
4089c2f3473SPeter Brune   ierr = PetscFree(jac->g);CHKERRQ(ierr);
40990a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
41090a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
41190a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
412eed5f15bSPeter Brune   PetscFunctionReturn(0);
413eed5f15bSPeter Brune }
414eed5f15bSPeter Brune 
415eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
416eed5f15bSPeter Brune {
417eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
418eed5f15bSPeter Brune   PetscErrorCode     ierr;
419eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
420eed5f15bSPeter Brune 
421eed5f15bSPeter Brune   PetscFunctionBegin;
422eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
423eed5f15bSPeter Brune   while (next) {
424eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
425eed5f15bSPeter Brune     next_tmp = next;
426eed5f15bSPeter Brune     next     = next->next;
427eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
428eed5f15bSPeter Brune   }
429eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
430eed5f15bSPeter Brune   PetscFunctionReturn(0);
431eed5f15bSPeter Brune }
432eed5f15bSPeter Brune 
4334416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,SNES snes)
434eed5f15bSPeter Brune {
435eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
436eed5f15bSPeter Brune   PetscErrorCode     ierr;
437eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
438eed5f15bSPeter Brune   SNES_CompositeLink next;
439eed5f15bSPeter Brune   char               *sneses[8];
4408f626970SPeter Brune   PetscReal          dmps[8];
441eed5f15bSPeter Brune   PetscBool          flg;
442eed5f15bSPeter Brune 
443eed5f15bSPeter Brune   PetscFunctionBegin;
444e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");CHKERRQ(ierr);
445eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
446eed5f15bSPeter Brune   if (flg) {
447eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
448eed5f15bSPeter Brune   }
449eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
450eed5f15bSPeter Brune   if (flg) {
451eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
452eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
453eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
454eed5f15bSPeter Brune     }
455eed5f15bSPeter Brune   }
4568f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
4578f626970SPeter Brune   if (flg) {
4588f626970SPeter Brune     for (i=0; i<nmax; i++) {
4598f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
4608f626970SPeter Brune     }
4618f626970SPeter Brune   }
4625e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr);
4635e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr);
464eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
465eed5f15bSPeter Brune 
466eed5f15bSPeter Brune   next = jac->head;
467eed5f15bSPeter Brune   while (next) {
468eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
469eed5f15bSPeter Brune     next = next->next;
470eed5f15bSPeter Brune   }
471eed5f15bSPeter Brune   PetscFunctionReturn(0);
472eed5f15bSPeter Brune }
473eed5f15bSPeter Brune 
474eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
475eed5f15bSPeter Brune {
476eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
477eed5f15bSPeter Brune   PetscErrorCode     ierr;
478eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
479eed5f15bSPeter Brune   PetscBool          iascii;
480eed5f15bSPeter Brune 
481eed5f15bSPeter Brune   PetscFunctionBegin;
482eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
483eed5f15bSPeter Brune   if (iascii) {
484efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
485eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
486eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n");CHKERRQ(ierr);
487eed5f15bSPeter Brune   }
488eed5f15bSPeter Brune   if (iascii) {
489eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
490eed5f15bSPeter Brune   }
491eed5f15bSPeter Brune   while (next) {
492eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
493eed5f15bSPeter Brune     next = next->next;
494eed5f15bSPeter Brune   }
495eed5f15bSPeter Brune   if (iascii) {
496eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
497eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n");CHKERRQ(ierr);
498eed5f15bSPeter Brune   }
499eed5f15bSPeter Brune   PetscFunctionReturn(0);
500eed5f15bSPeter Brune }
501eed5f15bSPeter Brune 
502eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
503eed5f15bSPeter Brune 
504eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
505eed5f15bSPeter Brune {
506eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
507eed5f15bSPeter Brune 
508eed5f15bSPeter Brune   PetscFunctionBegin;
509eed5f15bSPeter Brune   jac->type = type;
510eed5f15bSPeter Brune   PetscFunctionReturn(0);
511eed5f15bSPeter Brune }
512eed5f15bSPeter Brune 
513eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
514eed5f15bSPeter Brune {
515eed5f15bSPeter Brune   SNES_Composite     *jac;
516eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
517eed5f15bSPeter Brune   PetscErrorCode     ierr;
518eed5f15bSPeter Brune   PetscInt           cnt = 0;
519eed5f15bSPeter Brune   const char         *prefix;
520d726e3a5SJed Brown   char               newprefix[20];
521eed5f15bSPeter Brune   DM                 dm;
522eed5f15bSPeter Brune 
523eed5f15bSPeter Brune   PetscFunctionBegin;
524b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ilink);CHKERRQ(ierr);
525eed5f15bSPeter Brune   ilink->next = 0;
526eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
527eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
528eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
529eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
530cf5b3eb5SPeter Brune   ierr        = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr);
5315a94abe6SSatish Balay   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
532eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
533eed5f15bSPeter Brune   next = jac->head;
534eed5f15bSPeter Brune   if (!next) {
535eed5f15bSPeter Brune     jac->head       = ilink;
536eed5f15bSPeter Brune     ilink->previous = NULL;
537eed5f15bSPeter Brune   } else {
538eed5f15bSPeter Brune     cnt++;
539eed5f15bSPeter Brune     while (next->next) {
540eed5f15bSPeter Brune       next = next->next;
541eed5f15bSPeter Brune       cnt++;
542eed5f15bSPeter Brune     }
543eed5f15bSPeter Brune     next->next      = ilink;
544eed5f15bSPeter Brune     ilink->previous = next;
545eed5f15bSPeter Brune   }
546eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
547eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
548eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
549eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
5508f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
551eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
55272edecb9SPeter Brune   ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
55363cdc2bdSPatrick Farrell 
5548f626970SPeter Brune   ilink->dmp = 1.0;
55590a8ba9bSPeter Brune   jac->nsnes++;
556eed5f15bSPeter Brune   PetscFunctionReturn(0);
557eed5f15bSPeter Brune }
558eed5f15bSPeter Brune 
559eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
560eed5f15bSPeter Brune {
561eed5f15bSPeter Brune   SNES_Composite     *jac;
562eed5f15bSPeter Brune   SNES_CompositeLink next;
563eed5f15bSPeter Brune   PetscInt           i;
564eed5f15bSPeter Brune 
565eed5f15bSPeter Brune   PetscFunctionBegin;
566eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
567eed5f15bSPeter Brune   next = jac->head;
568eed5f15bSPeter Brune   for (i=0; i<n; i++) {
569eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
570eed5f15bSPeter Brune     next = next->next;
571eed5f15bSPeter Brune   }
572eed5f15bSPeter Brune   *subsnes = next->snes;
573eed5f15bSPeter Brune   PetscFunctionReturn(0);
574eed5f15bSPeter Brune }
575eed5f15bSPeter Brune 
576eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
577e31ab4f9SPeter Brune /*@C
578eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
579eed5f15bSPeter Brune 
580eed5f15bSPeter Brune    Logically Collective on SNES
581eed5f15bSPeter Brune 
582eed5f15bSPeter Brune    Input Parameter:
583eed5f15bSPeter Brune +  snes - the preconditioner context
584eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
585eed5f15bSPeter Brune 
586eed5f15bSPeter Brune    Options Database Key:
587eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
588eed5f15bSPeter Brune 
589eed5f15bSPeter Brune    Level: Developer
590eed5f15bSPeter Brune 
591eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
592eed5f15bSPeter Brune @*/
593eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
594eed5f15bSPeter Brune {
595eed5f15bSPeter Brune   PetscErrorCode ierr;
596eed5f15bSPeter Brune 
597eed5f15bSPeter Brune   PetscFunctionBegin;
598eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
599eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
600eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
601eed5f15bSPeter Brune   PetscFunctionReturn(0);
602eed5f15bSPeter Brune }
603eed5f15bSPeter Brune 
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 /*@
627eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
628eed5f15bSPeter Brune 
629eed5f15bSPeter Brune    Not Collective
630eed5f15bSPeter Brune 
631eed5f15bSPeter Brune    Input Parameter:
632eed5f15bSPeter Brune +  snes - the preconditioner context
633eed5f15bSPeter Brune -  n - the number of the snes requested
634eed5f15bSPeter Brune 
635eed5f15bSPeter Brune    Output Parameters:
636eed5f15bSPeter Brune .  subsnes - the SNES requested
637eed5f15bSPeter Brune 
638eed5f15bSPeter Brune    Level: Developer
639eed5f15bSPeter Brune 
640eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
641eed5f15bSPeter Brune 
642eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
643eed5f15bSPeter Brune @*/
644eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
645eed5f15bSPeter Brune {
646eed5f15bSPeter Brune   PetscErrorCode ierr;
647eed5f15bSPeter Brune 
648eed5f15bSPeter Brune   PetscFunctionBegin;
649eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
650eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
651eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
652eed5f15bSPeter Brune   PetscFunctionReturn(0);
653eed5f15bSPeter Brune }
654eed5f15bSPeter Brune 
6556b2b7f7bSPatrick Farrell /*@
6566b2b7f7bSPatrick Farrell    SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES.
6576b2b7f7bSPatrick Farrell 
6586b2b7f7bSPatrick Farrell    Logically Collective on SNES
6596b2b7f7bSPatrick Farrell 
6606b2b7f7bSPatrick Farrell    Input Parameter:
6616b2b7f7bSPatrick Farrell    snes - the preconditioner context
6626b2b7f7bSPatrick Farrell 
6636b2b7f7bSPatrick Farrell    Output Parameter:
6646b2b7f7bSPatrick Farrell    n - the number of subsolvers
6656b2b7f7bSPatrick Farrell 
6666b2b7f7bSPatrick Farrell    Level: Developer
6676b2b7f7bSPatrick Farrell 
6686b2b7f7bSPatrick Farrell .keywords: SNES, composite preconditioner
6696b2b7f7bSPatrick Farrell @*/
6706b2b7f7bSPatrick Farrell PetscErrorCode  SNESCompositeGetNumber(SNES snes,PetscInt *n)
6716b2b7f7bSPatrick Farrell {
6726b2b7f7bSPatrick Farrell   SNES_Composite     *jac;
6736b2b7f7bSPatrick Farrell   SNES_CompositeLink next;
6746b2b7f7bSPatrick Farrell 
6756b2b7f7bSPatrick Farrell   PetscFunctionBegin;
6766b2b7f7bSPatrick Farrell   jac  = (SNES_Composite*)snes->data;
6776b2b7f7bSPatrick Farrell   next = jac->head;
6786b2b7f7bSPatrick Farrell 
6796b2b7f7bSPatrick Farrell   *n = 0;
6806b2b7f7bSPatrick Farrell   while (next) {
6816b2b7f7bSPatrick Farrell     *n = *n + 1;
6826b2b7f7bSPatrick Farrell     next = next->next;
6836b2b7f7bSPatrick Farrell   }
6846b2b7f7bSPatrick Farrell   PetscFunctionReturn(0);
6856b2b7f7bSPatrick Farrell }
6866b2b7f7bSPatrick Farrell 
6878f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
6888f626970SPeter Brune {
6898f626970SPeter Brune   SNES_Composite     *jac;
6908f626970SPeter Brune   SNES_CompositeLink next;
6918f626970SPeter Brune   PetscInt           i;
6928f626970SPeter Brune 
6938f626970SPeter Brune   PetscFunctionBegin;
6948f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
6958f626970SPeter Brune   next = jac->head;
6968f626970SPeter Brune   for (i=0; i<n; i++) {
6978f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
6988f626970SPeter Brune     next = next->next;
6998f626970SPeter Brune   }
7008f626970SPeter Brune   next->dmp = dmp;
7018f626970SPeter Brune   PetscFunctionReturn(0);
7028f626970SPeter Brune }
7038f626970SPeter Brune 
7048f626970SPeter Brune /*@
7058f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
7068f626970SPeter Brune 
7078f626970SPeter Brune    Not Collective
7088f626970SPeter Brune 
7098f626970SPeter Brune    Input Parameter:
7108f626970SPeter Brune +  snes - the preconditioner context
7118f626970SPeter Brune .  n - the number of the snes requested
7128f626970SPeter Brune -  dmp - the damping
7138f626970SPeter Brune 
7148f626970SPeter Brune    Level: Developer
7158f626970SPeter Brune 
7168f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
7178f626970SPeter Brune 
7188f626970SPeter Brune .seealso: SNESCompositeAddSNES()
7198f626970SPeter Brune @*/
7208f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
7218f626970SPeter Brune {
7228f626970SPeter Brune   PetscErrorCode ierr;
7238f626970SPeter Brune 
7248f626970SPeter Brune   PetscFunctionBegin;
7258f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7268f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
7278f626970SPeter Brune   PetscFunctionReturn(0);
7288f626970SPeter Brune }
7298f626970SPeter Brune 
73040244768SBarry Smith static PetscErrorCode SNESSolve_Composite(SNES snes)
731eed5f15bSPeter Brune {
732eed5f15bSPeter Brune   Vec              F;
733eed5f15bSPeter Brune   Vec              X;
734eed5f15bSPeter Brune   Vec              B;
735eed5f15bSPeter Brune   PetscInt         i;
736b22975d2SPatrick Farrell   PetscReal        fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
737eed5f15bSPeter Brune   PetscErrorCode   ierr;
73872edecb9SPeter Brune   SNESNormSchedule normtype;
739eed5f15bSPeter Brune   SNES_Composite   *comp = (SNES_Composite*)snes->data;
740eed5f15bSPeter Brune 
741eed5f15bSPeter Brune   PetscFunctionBegin;
742eed5f15bSPeter Brune   X = snes->vec_sol;
743eed5f15bSPeter Brune   F = snes->vec_func;
744eed5f15bSPeter Brune   B = snes->vec_rhs;
745eed5f15bSPeter Brune 
746e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
747eed5f15bSPeter Brune   snes->iter   = 0;
748eed5f15bSPeter Brune   snes->norm   = 0.;
7490b762f1fSPatrick Farrell   comp->innerFailures = 0;
750e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
751b22975d2SPatrick Farrell   ierr         = SNESSetWorkVecs(snes, 1);CHKERRQ(ierr);
752eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
75372edecb9SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr);
75472edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
755eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
756eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
757eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
758eed5f15bSPeter Brune 
759a6da83ecSPatrick Farrell     if (snes->xl && snes->xu) {
760a6da83ecSPatrick Farrell       ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr);
761a6da83ecSPatrick Farrell     } else {
762eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
763a6da83ecSPatrick Farrell     }
764422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
765e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
766eed5f15bSPeter Brune     snes->iter = 0;
767eed5f15bSPeter Brune     snes->norm = fnorm;
768e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
769eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
770eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
771eed5f15bSPeter Brune 
772eed5f15bSPeter Brune     /* test convergence */
773eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
774eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
775eed5f15bSPeter Brune   } else {
776e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
777eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
778eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
779eed5f15bSPeter Brune   }
780eed5f15bSPeter Brune 
781eed5f15bSPeter Brune   /* Call general purpose update function */
782eed5f15bSPeter Brune   if (snes->ops->update) {
783eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
784eed5f15bSPeter Brune   }
785eed5f15bSPeter Brune 
786eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
787b22975d2SPatrick Farrell     /* Copy the state before modification by application of the composite solver;
788b22975d2SPatrick Farrell        we will subtract the new state after application */
789b22975d2SPatrick Farrell     ierr = VecCopy(X, snes->work[0]);CHKERRQ(ierr);
790b22975d2SPatrick Farrell 
791eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
792eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
793eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
794eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
79590a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
79690a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
7976c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
798d5a53f06SPeter Brune     if (snes->reason < 0) break;
799d5a53f06SPeter Brune 
800b22975d2SPatrick Farrell     /* Compute the solution update for convergence testing */
801b22975d2SPatrick Farrell     ierr = VecAXPY(snes->work[0], -1.0, X);CHKERRQ(ierr);
802b22975d2SPatrick Farrell     ierr = VecScale(snes->work[0], -1.0);CHKERRQ(ierr);
803b22975d2SPatrick Farrell 
804eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
805eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
806b22975d2SPatrick Farrell 
807a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
808a6da83ecSPatrick Farrell         ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
809a6da83ecSPatrick Farrell         ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
810a6da83ecSPatrick Farrell         ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr);
811a6da83ecSPatrick Farrell         ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
812a6da83ecSPatrick Farrell         ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
813a6da83ecSPatrick Farrell       } else {
814b22975d2SPatrick Farrell         ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr);
815b22975d2SPatrick Farrell         ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
816b22975d2SPatrick Farrell         ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
817b22975d2SPatrick Farrell 
818b22975d2SPatrick Farrell         ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr);
819b22975d2SPatrick Farrell         ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
820b22975d2SPatrick Farrell         ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
821a6da83ecSPatrick Farrell       }
822422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
823b22975d2SPatrick Farrell     } else if (normtype == SNES_NORM_ALWAYS) {
824b22975d2SPatrick Farrell       ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
825b22975d2SPatrick Farrell       ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
826b22975d2SPatrick Farrell       ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
827b22975d2SPatrick Farrell       ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr);
828eed5f15bSPeter Brune     }
829eed5f15bSPeter Brune     /* Monitor convergence */
830e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
831eed5f15bSPeter Brune     snes->iter = i+1;
832eed5f15bSPeter Brune     snes->norm = fnorm;
833*c1e67a49SFande Kong     snes->xnorm = xnorm;
834*c1e67a49SFande Kong     snes->ynorm = snorm;
835e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
836eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
837eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
838eed5f15bSPeter Brune     /* Test for convergence */
839b22975d2SPatrick Farrell     if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
840eed5f15bSPeter Brune     if (snes->reason) break;
841eed5f15bSPeter Brune     /* Call general purpose update function */
842eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
843eed5f15bSPeter Brune   }
84472edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS) {
845eed5f15bSPeter Brune     if (i == snes->max_its) {
846eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
847eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
848eed5f15bSPeter Brune     }
849eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
850eed5f15bSPeter Brune   PetscFunctionReturn(0);
851eed5f15bSPeter Brune }
852eed5f15bSPeter Brune 
853eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
854eed5f15bSPeter Brune 
855eed5f15bSPeter Brune /*MC
856eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
857eed5f15bSPeter Brune 
858eed5f15bSPeter Brune    Options Database Keys:
859eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
860eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
861eed5f15bSPeter Brune 
862eed5f15bSPeter Brune    Level: intermediate
863eed5f15bSPeter Brune 
864eed5f15bSPeter Brune    Concepts: composing solvers
865eed5f15bSPeter Brune 
866eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
867eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
868eed5f15bSPeter Brune            SNESCompositeGetSNES()
869eed5f15bSPeter Brune 
8704f02bc6aSBarry Smith    References:
87196a0c994SBarry Smith .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
8724f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
8734f02bc6aSBarry Smith 
874eed5f15bSPeter Brune M*/
875eed5f15bSPeter Brune 
876eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
877eed5f15bSPeter Brune {
878eed5f15bSPeter Brune   PetscErrorCode ierr;
879eed5f15bSPeter Brune   SNES_Composite *jac;
880eed5f15bSPeter Brune 
881eed5f15bSPeter Brune   PetscFunctionBegin;
882b00a9115SJed Brown   ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr);
883eed5f15bSPeter Brune 
884eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
885eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
886eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
887eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
888eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
889eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
890eed5f15bSPeter Brune 
8914fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8924fc747eaSLawrence Mitchell 
893eed5f15bSPeter Brune   snes->data = (void*)jac;
89490a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
89590a8ba9bSPeter Brune   jac->Fes   = NULL;
89690a8ba9bSPeter Brune   jac->Xes   = NULL;
8979c2f3473SPeter Brune   jac->fnorms = NULL;
89890a8ba9bSPeter Brune   jac->nsnes = 0;
899eed5f15bSPeter Brune   jac->head  = 0;
9005e806d2eSPeter Brune   jac->stol  = 0.1;
9015e806d2eSPeter Brune   jac->rtol  = 1.1;
902eed5f15bSPeter Brune 
90390a8ba9bSPeter Brune   jac->h     = NULL;
90490a8ba9bSPeter Brune   jac->s     = NULL;
90590a8ba9bSPeter Brune   jac->beta  = NULL;
90690a8ba9bSPeter Brune   jac->work  = NULL;
90790a8ba9bSPeter Brune   jac->rwork = NULL;
90890a8ba9bSPeter Brune 
9098f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
9108f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
9118f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
9128f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
913eed5f15bSPeter Brune   PetscFunctionReturn(0);
914eed5f15bSPeter Brune }
915