xref: /petsc/src/snes/impls/composite/snescomposite.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
89e5d0892SLisandro Dalcin const char *const        SNESCompositeTypes[]   = {"ADDITIVE","MULTIPLICATIVE","ADDITIVEOPTIMAL","SNESCompositeType","SNES_COMPOSITE",NULL};
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   SNES_Composite      *jac = (SNES_Composite*)snes->data;
51eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
52eed5f15bSPeter Brune   Vec                 FSub;
53d5a53f06SPeter Brune   SNESConvergedReason reason;
54eed5f15bSPeter Brune 
55eed5f15bSPeter Brune   PetscFunctionBegin;
5628b400f6SJacob Faibussowitsch   PetscCheck(next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
5772edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
589566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next->snes,F));
59eed5f15bSPeter Brune   }
609566063dSJacob Faibussowitsch   PetscCall(SNESSolve(next->snes,B,X));
619566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(next->snes,&reason));
62d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
630b762f1fSPatrick Farrell     jac->innerFailures++;
640b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
65d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
66d5a53f06SPeter Brune       PetscFunctionReturn(0);
67d5a53f06SPeter Brune     }
680b762f1fSPatrick Farrell   }
69eed5f15bSPeter Brune 
70eed5f15bSPeter Brune   while (next->next) {
71eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
72efd4aadfSBarry Smith     if (next->snes->npcside== PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
739566063dSJacob Faibussowitsch       PetscCall(SNESGetFunction(next->snes,&FSub,NULL,NULL));
74eed5f15bSPeter Brune       next = next->next;
759566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(next->snes,FSub));
76eed5f15bSPeter Brune     } else {
77eed5f15bSPeter Brune       next = next->next;
78eed5f15bSPeter Brune     }
799566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next->snes,B,X));
809566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next->snes,&reason));
81d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
820b762f1fSPatrick Farrell       jac->innerFailures++;
830b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
84d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
85d5a53f06SPeter Brune         PetscFunctionReturn(0);
86d5a53f06SPeter Brune       }
87eed5f15bSPeter Brune     }
880b762f1fSPatrick Farrell   }
89efd4aadfSBarry Smith   if (next->snes->npcside== PC_RIGHT) {
909566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(next->snes,&FSub,NULL,NULL));
919566063dSJacob Faibussowitsch     PetscCall(VecCopy(FSub,F));
92d5a53f06SPeter Brune     if (fnorm) {
9363cdc2bdSPatrick Farrell       if (snes->xl && snes->xu) {
949566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
9563cdc2bdSPatrick Farrell       } else {
969566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, fnorm));
9763cdc2bdSPatrick Farrell       }
98422a814eSBarry Smith       SNESCheckFunctionNorm(snes,*fnorm);
99d5a53f06SPeter Brune     }
10072edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
1019566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes,X,F));
102d5a53f06SPeter Brune     if (fnorm) {
103a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
1049566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
105a6da83ecSPatrick Farrell       } else {
1069566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, fnorm));
107a6da83ecSPatrick Farrell       }
108422a814eSBarry Smith       SNESCheckFunctionNorm(snes,*fnorm);
109d5a53f06SPeter Brune     }
110eed5f15bSPeter Brune   }
111eed5f15bSPeter Brune   PetscFunctionReturn(0);
112eed5f15bSPeter Brune }
113eed5f15bSPeter Brune 
114eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
115eed5f15bSPeter Brune {
116eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
117eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
118eed5f15bSPeter Brune   Vec                 Y,Xorig;
119d5a53f06SPeter Brune   SNESConvergedReason reason;
120eed5f15bSPeter Brune 
121eed5f15bSPeter Brune   PetscFunctionBegin;
122eed5f15bSPeter Brune   Y = snes->vec_sol_update;
1239566063dSJacob Faibussowitsch   if (!jac->Xorig) PetscCall(VecDuplicate(X,&jac->Xorig));
124eed5f15bSPeter Brune   Xorig = jac->Xorig;
1259566063dSJacob Faibussowitsch   PetscCall(VecCopy(X,Xorig));
12628b400f6SJacob Faibussowitsch   PetscCheck(next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
12772edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
1289566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next->snes,F));
129eed5f15bSPeter Brune     while (next->next) {
130eed5f15bSPeter Brune       next = next->next;
1319566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(next->snes,F));
132eed5f15bSPeter Brune     }
133eed5f15bSPeter Brune   }
134eed5f15bSPeter Brune   next = jac->head;
1359566063dSJacob Faibussowitsch   PetscCall(VecCopy(Xorig,Y));
1369566063dSJacob Faibussowitsch   PetscCall(SNESSolve(next->snes,B,Y));
1379566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(next->snes,&reason));
138d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1390b762f1fSPatrick Farrell     jac->innerFailures++;
1400b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
141d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
142d5a53f06SPeter Brune       PetscFunctionReturn(0);
143d5a53f06SPeter Brune     }
1440b762f1fSPatrick Farrell   }
1459566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y,-1.0,Xorig));
1469566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,next->dmp,Y));
1478f626970SPeter Brune   while (next->next) {
1488f626970SPeter Brune     next = next->next;
1499566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xorig,Y));
1509566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next->snes,B,Y));
1519566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next->snes,&reason));
152d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1530b762f1fSPatrick Farrell       jac->innerFailures++;
1540b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
155d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
156d5a53f06SPeter Brune         PetscFunctionReturn(0);
157d5a53f06SPeter Brune       }
1580b762f1fSPatrick Farrell     }
1599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y,-1.0,Xorig));
1609566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X,next->dmp,Y));
161eed5f15bSPeter Brune   }
16272edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
1639566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes,X,F));
164a6da83ecSPatrick Farrell     if (fnorm) {
165a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
1669566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
167a6da83ecSPatrick Farrell       } else {
1689566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, fnorm));
169a6da83ecSPatrick Farrell       }
170422a814eSBarry Smith       SNESCheckFunctionNorm(snes,*fnorm);
171a6da83ecSPatrick Farrell     }
172eed5f15bSPeter Brune   }
173eed5f15bSPeter Brune   PetscFunctionReturn(0);
174eed5f15bSPeter Brune }
17590a8ba9bSPeter Brune 
17690a8ba9bSPeter Brune /* approximately solve the overdetermined system:
17790a8ba9bSPeter Brune 
17890a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
17990a8ba9bSPeter Brune  \alpha_i                      = 1
18090a8ba9bSPeter Brune 
18190a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
18290a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
18390a8ba9bSPeter Brune 
18490a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
18590a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
18690a8ba9bSPeter Brune  */
18790a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
18890a8ba9bSPeter Brune {
18990a8ba9bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
19090a8ba9bSPeter Brune   SNES_CompositeLink  next = jac->head;
19190a8ba9bSPeter Brune   Vec                 *Xes = jac->Xes,*Fes = jac->Fes;
19290a8ba9bSPeter Brune   PetscInt            i,j;
1935e806d2eSPeter Brune   PetscScalar         tot,total,ftf;
1945e806d2eSPeter Brune   PetscReal           min_fnorm;
1955e806d2eSPeter Brune   PetscInt            min_i;
196d5a53f06SPeter Brune   SNESConvergedReason reason;
19790a8ba9bSPeter Brune 
19890a8ba9bSPeter Brune   PetscFunctionBegin;
19928b400f6SJacob Faibussowitsch   PetscCheck(next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
20058bdfa14SPeter Brune 
20158bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
20258bdfa14SPeter Brune     next = jac->head;
2039566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next->snes,F));
20458bdfa14SPeter Brune     while (next->next) {
20558bdfa14SPeter Brune       next = next->next;
2069566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(next->snes,F));
20758bdfa14SPeter Brune     }
20858bdfa14SPeter Brune   }
20958bdfa14SPeter Brune 
21090a8ba9bSPeter Brune   next = jac->head;
21190a8ba9bSPeter Brune   i = 0;
2129566063dSJacob Faibussowitsch   PetscCall(VecCopy(X,Xes[i]));
2139566063dSJacob Faibussowitsch   PetscCall(SNESSolve(next->snes,B,Xes[i]));
2149566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(next->snes,&reason));
215d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2160b762f1fSPatrick Farrell     jac->innerFailures++;
2170b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
218d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
219d5a53f06SPeter Brune       PetscFunctionReturn(0);
220d5a53f06SPeter Brune     }
2210b762f1fSPatrick Farrell   }
22290a8ba9bSPeter Brune   while (next->next) {
22390a8ba9bSPeter Brune     i++;
22490a8ba9bSPeter Brune     next = next->next;
2259566063dSJacob Faibussowitsch     PetscCall(VecCopy(X,Xes[i]));
2269566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next->snes,B,Xes[i]));
2279566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next->snes,&reason));
228d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2290b762f1fSPatrick Farrell       jac->innerFailures++;
2300b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
231d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
232d5a53f06SPeter Brune         PetscFunctionReturn(0);
233d5a53f06SPeter Brune       }
23490a8ba9bSPeter Brune     }
2350b762f1fSPatrick Farrell   }
23690a8ba9bSPeter Brune 
23790a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
23890a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
23990a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2409566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]));
24190a8ba9bSPeter Brune     }
2429566063dSJacob Faibussowitsch     PetscCall(VecDotBegin(Fes[i],F,&jac->g[i]));
24390a8ba9bSPeter Brune   }
2445e806d2eSPeter Brune 
24590a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
24690a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2479566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]));
2489c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
24990a8ba9bSPeter Brune     }
2509566063dSJacob Faibussowitsch     PetscCall(VecDotEnd(Fes[i],F,&jac->g[i]));
25190a8ba9bSPeter Brune   }
25290a8ba9bSPeter Brune 
2539c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
25490a8ba9bSPeter Brune 
25590a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
25690a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
2579c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
25890a8ba9bSPeter Brune     }
25990a8ba9bSPeter Brune   }
26090a8ba9bSPeter Brune 
26190a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
2629c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
2639c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
26490a8ba9bSPeter Brune     }
2659c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2669c2f3473SPeter Brune   }
26790a8ba9bSPeter Brune 
26890a8ba9bSPeter Brune   jac->info  = 0;
26990a8ba9bSPeter Brune   jac->rcond = -1.;
2709566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
27190a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
272792fecdfSBarry Smith   PetscCallBLAS("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));
27390a8ba9bSPeter Brune #else
274792fecdfSBarry Smith   PetscCallBLAS("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));
27590a8ba9bSPeter Brune #endif
2769566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
27708401ef6SPierre Jolivet   PetscCheck(jac->info >= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
27808401ef6SPierre Jolivet   PetscCheck(jac->info <= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
2799c2f3473SPeter Brune   tot = 0.;
2805e806d2eSPeter Brune   total = 0.;
28190a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
28208401ef6SPierre Jolivet     PetscCheck(!snes->errorifnotconverged || !PetscIsInfOrNanScalar(jac->beta[i]),PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
28363a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"%" PetscInt_FMT ": %g\n",i,(double)PetscRealPart(jac->beta[i])));
2849c2f3473SPeter Brune     tot += jac->beta[i];
2855e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
28690a8ba9bSPeter Brune   }
2879566063dSJacob Faibussowitsch   PetscCall(VecScale(X,(1. - tot)));
2889566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X,jac->n,jac->beta,Xes));
2899566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes,X,F));
290a6da83ecSPatrick Farrell 
291a6da83ecSPatrick Farrell   if (snes->xl && snes->xu) {
2929566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
293a6da83ecSPatrick Farrell   } else {
2949566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, fnorm));
295a6da83ecSPatrick Farrell   }
29690a8ba9bSPeter Brune 
2975e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
2985e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
2995e806d2eSPeter Brune   min_i     = 0;
3009c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
3015e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
3025e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
3035e806d2eSPeter Brune       min_i     = i;
3049c2f3473SPeter Brune     }
3059c2f3473SPeter Brune   }
3065e806d2eSPeter Brune 
3075e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
3085e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
3099566063dSJacob Faibussowitsch     PetscCall(VecCopy(jac->Xes[min_i],X));
3109566063dSJacob Faibussowitsch     PetscCall(VecCopy(jac->Fes[min_i],F));
3115e806d2eSPeter Brune     *fnorm = min_fnorm;
3125e806d2eSPeter Brune   }
31390a8ba9bSPeter Brune   PetscFunctionReturn(0);
31490a8ba9bSPeter Brune }
31590a8ba9bSPeter Brune 
316eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
317eed5f15bSPeter Brune {
318dd515a93SPeter Brune   DM                 dm;
319eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
320eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
32190a8ba9bSPeter Brune   PetscInt           n=0,i;
32290a8ba9bSPeter Brune   Vec                F;
323eed5f15bSPeter Brune 
324eed5f15bSPeter Brune   PetscFunctionBegin;
3259566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
32663cdc2bdSPatrick Farrell 
32763cdc2bdSPatrick Farrell   if (snes->ops->computevariablebounds) {
32863cdc2bdSPatrick Farrell     /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */
3299566063dSJacob Faibussowitsch     if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol,&snes->xl));
3309566063dSJacob Faibussowitsch     if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol,&snes->xu));
331*dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes,computevariablebounds ,snes->xl,snes->xu);
33263cdc2bdSPatrick Farrell   }
33363cdc2bdSPatrick Farrell 
334eed5f15bSPeter Brune   while (next) {
33590a8ba9bSPeter Brune     n++;
3369566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(next->snes,dm));
3379566063dSJacob Faibussowitsch     PetscCall(SNESSetApplicationContext(next->snes, snes->user));
33863cdc2bdSPatrick Farrell     if (snes->xl && snes->xu) {
33963cdc2bdSPatrick Farrell       if (snes->ops->computevariablebounds) {
3409566063dSJacob Faibussowitsch         PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds));
34163cdc2bdSPatrick Farrell       } else {
3429566063dSJacob Faibussowitsch         PetscCall(SNESVISetVariableBounds(next->snes,snes->xl,snes->xu));
34363cdc2bdSPatrick Farrell       }
34463cdc2bdSPatrick Farrell     }
34563cdc2bdSPatrick Farrell 
346eed5f15bSPeter Brune     next = next->next;
347eed5f15bSPeter Brune   }
34890a8ba9bSPeter Brune   jac->nsnes = n;
3499566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes,&F,NULL,NULL));
35090a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
3519566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(F,jac->nsnes,&jac->Xes));
3529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&jac->Fes));
3539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&jac->fnorms));
35490a8ba9bSPeter Brune     next = jac->head;
35590a8ba9bSPeter Brune     i = 0;
35690a8ba9bSPeter Brune     while (next) {
3579566063dSJacob Faibussowitsch       PetscCall(SNESGetFunction(next->snes,&F,NULL,NULL));
35890a8ba9bSPeter Brune       jac->Fes[i] = F;
3599566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)F));
36090a8ba9bSPeter Brune       next = next->next;
36190a8ba9bSPeter Brune       i++;
36290a8ba9bSPeter Brune     }
36390a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
36490a8ba9bSPeter Brune     jac->nrhs  = 1;
3659c2f3473SPeter Brune     jac->lda   = jac->nsnes;
36690a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
36790a8ba9bSPeter Brune     jac->n     = jac->nsnes;
36890a8ba9bSPeter Brune 
3699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n*jac->n,&jac->h));
3709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n,&jac->beta));
3719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n,&jac->s));
3729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n,&jac->g));
3739c2f3473SPeter Brune     jac->lwork = 12*jac->n;
374dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
3759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->lwork,&jac->rwork));
37690a8ba9bSPeter Brune #endif
3779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->lwork,&jac->work));
37890a8ba9bSPeter Brune   }
37990a8ba9bSPeter Brune 
380eed5f15bSPeter Brune   PetscFunctionReturn(0);
381eed5f15bSPeter Brune }
382eed5f15bSPeter Brune 
383eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
384eed5f15bSPeter Brune {
385eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
386eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
387eed5f15bSPeter Brune 
388eed5f15bSPeter Brune   PetscFunctionBegin;
389eed5f15bSPeter Brune   while (next) {
3909566063dSJacob Faibussowitsch     PetscCall(SNESReset(next->snes));
391eed5f15bSPeter Brune     next = next->next;
392eed5f15bSPeter Brune   }
3939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->Xorig));
3949566063dSJacob Faibussowitsch   if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes,&jac->Xes));
3959566063dSJacob Faibussowitsch   if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes,&jac->Fes));
3969566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->fnorms));
3979566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->h));
3989566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->s));
3999566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g));
4009566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->beta));
4019566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->work));
4029566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->rwork));
403eed5f15bSPeter Brune   PetscFunctionReturn(0);
404eed5f15bSPeter Brune }
405eed5f15bSPeter Brune 
406eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
407eed5f15bSPeter Brune {
408eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
409eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
410eed5f15bSPeter Brune 
411eed5f15bSPeter Brune   PetscFunctionBegin;
4129566063dSJacob Faibussowitsch   PetscCall(SNESReset_Composite(snes));
413eed5f15bSPeter Brune   while (next) {
4149566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&next->snes));
415eed5f15bSPeter Brune     next_tmp = next;
416eed5f15bSPeter Brune     next     = next->next;
4179566063dSJacob Faibussowitsch     PetscCall(PetscFree(next_tmp));
418eed5f15bSPeter Brune   }
4192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",NULL));
4202e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",NULL));
4212e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",NULL));
4222e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",NULL));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
424eed5f15bSPeter Brune   PetscFunctionReturn(0);
425eed5f15bSPeter Brune }
426eed5f15bSPeter Brune 
427*dbbe0bcdSBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(SNES snes,PetscOptionItems *PetscOptionsObject)
428eed5f15bSPeter Brune {
429eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
430eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
431eed5f15bSPeter Brune   SNES_CompositeLink next;
432eed5f15bSPeter Brune   char               *sneses[8];
4338f626970SPeter Brune   PetscReal          dmps[8];
434eed5f15bSPeter Brune   PetscBool          flg;
435eed5f15bSPeter Brune 
436eed5f15bSPeter Brune   PetscFunctionBegin;
437d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Composite preconditioner options");
4389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg));
4391baa6e33SBarry Smith   if (flg) PetscCall(SNESCompositeSetType(snes,jac->type));
4409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg));
441eed5f15bSPeter Brune   if (flg) {
442eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
4439566063dSJacob Faibussowitsch       PetscCall(SNESCompositeAddSNES(snes,sneses[i]));
4449566063dSJacob Faibussowitsch       PetscCall(PetscFree(sneses[i]));   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
445eed5f15bSPeter Brune     }
446eed5f15bSPeter Brune   }
4479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg));
4488f626970SPeter Brune   if (flg) {
4498f626970SPeter Brune     for (i=0; i<nmax; i++) {
4509566063dSJacob Faibussowitsch       PetscCall(SNESCompositeSetDamping(snes,i,dmps[i]));
4518f626970SPeter Brune     }
4528f626970SPeter Brune   }
4539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL));
4549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL));
455d0609cedSBarry Smith   PetscOptionsHeadEnd();
456eed5f15bSPeter Brune 
457eed5f15bSPeter Brune   next = jac->head;
458eed5f15bSPeter Brune   while (next) {
4599566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(next->snes));
460eed5f15bSPeter Brune     next = next->next;
461eed5f15bSPeter Brune   }
462eed5f15bSPeter Brune   PetscFunctionReturn(0);
463eed5f15bSPeter Brune }
464eed5f15bSPeter Brune 
465eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
466eed5f15bSPeter Brune {
467eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
468eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
469eed5f15bSPeter Brune   PetscBool          iascii;
470eed5f15bSPeter Brune 
471eed5f15bSPeter Brune   PetscFunctionBegin;
4729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
473eed5f15bSPeter Brune   if (iascii) {
4749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  type - %s\n",SNESCompositeTypes[jac->type]));
4759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  SNESes on composite preconditioner follow\n"));
4769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n"));
477eed5f15bSPeter Brune   }
4781baa6e33SBarry Smith   if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer));
479eed5f15bSPeter Brune   while (next) {
4809566063dSJacob Faibussowitsch     PetscCall(SNESView(next->snes,viewer));
481eed5f15bSPeter Brune     next = next->next;
482eed5f15bSPeter Brune   }
483eed5f15bSPeter Brune   if (iascii) {
4849566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
4859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n"));
486eed5f15bSPeter Brune   }
487eed5f15bSPeter Brune   PetscFunctionReturn(0);
488eed5f15bSPeter Brune }
489eed5f15bSPeter Brune 
490eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
491eed5f15bSPeter Brune 
492eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
493eed5f15bSPeter Brune {
494eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
495eed5f15bSPeter Brune 
496eed5f15bSPeter Brune   PetscFunctionBegin;
497eed5f15bSPeter Brune   jac->type = type;
498eed5f15bSPeter Brune   PetscFunctionReturn(0);
499eed5f15bSPeter Brune }
500eed5f15bSPeter Brune 
501eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
502eed5f15bSPeter Brune {
503eed5f15bSPeter Brune   SNES_Composite     *jac;
504eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
505eed5f15bSPeter Brune   PetscInt           cnt = 0;
506eed5f15bSPeter Brune   const char         *prefix;
507d726e3a5SJed Brown   char               newprefix[20];
508eed5f15bSPeter Brune   DM                 dm;
509eed5f15bSPeter Brune 
510eed5f15bSPeter Brune   PetscFunctionBegin;
5119566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&ilink));
5129e5d0892SLisandro Dalcin   ilink->next = NULL;
5139566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes));
5149566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes));
5159566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
5169566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(ilink->snes,dm));
5179566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(ilink->snes,snes->abstol,snes->rtol,snes->stol,1,snes->max_funcs));
5189566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)ilink->snes));
519eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
520eed5f15bSPeter Brune   next = jac->head;
521eed5f15bSPeter Brune   if (!next) {
522eed5f15bSPeter Brune     jac->head       = ilink;
523eed5f15bSPeter Brune     ilink->previous = NULL;
524eed5f15bSPeter Brune   } else {
525eed5f15bSPeter Brune     cnt++;
526eed5f15bSPeter Brune     while (next->next) {
527eed5f15bSPeter Brune       next = next->next;
528eed5f15bSPeter Brune       cnt++;
529eed5f15bSPeter Brune     }
530eed5f15bSPeter Brune     next->next      = ilink;
531eed5f15bSPeter Brune     ilink->previous = next;
532eed5f15bSPeter Brune   }
5339566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes,&prefix));
5349566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(ilink->snes,prefix));
5359566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(newprefix,sizeof(newprefix),"sub_%d_",(int)cnt));
5369566063dSJacob Faibussowitsch   PetscCall(SNESAppendOptionsPrefix(ilink->snes,newprefix));
5379566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1));
5389566063dSJacob Faibussowitsch   PetscCall(SNESSetType(ilink->snes,type));
5399566063dSJacob Faibussowitsch   PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY));
54063cdc2bdSPatrick Farrell 
5418f626970SPeter Brune   ilink->dmp = 1.0;
54290a8ba9bSPeter Brune   jac->nsnes++;
543eed5f15bSPeter Brune   PetscFunctionReturn(0);
544eed5f15bSPeter Brune }
545eed5f15bSPeter Brune 
546eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
547eed5f15bSPeter Brune {
548eed5f15bSPeter Brune   SNES_Composite     *jac;
549eed5f15bSPeter Brune   SNES_CompositeLink next;
550eed5f15bSPeter Brune   PetscInt           i;
551eed5f15bSPeter Brune 
552eed5f15bSPeter Brune   PetscFunctionBegin;
553eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
554eed5f15bSPeter Brune   next = jac->head;
555eed5f15bSPeter Brune   for (i=0; i<n; i++) {
55628b400f6SJacob Faibussowitsch     PetscCheck(next->next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
557eed5f15bSPeter Brune     next = next->next;
558eed5f15bSPeter Brune   }
559eed5f15bSPeter Brune   *subsnes = next->snes;
560eed5f15bSPeter Brune   PetscFunctionReturn(0);
561eed5f15bSPeter Brune }
562eed5f15bSPeter Brune 
563eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
564e31ab4f9SPeter Brune /*@C
565eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
566eed5f15bSPeter Brune 
567eed5f15bSPeter Brune    Logically Collective on SNES
568eed5f15bSPeter Brune 
569d8d19677SJose E. Roman    Input Parameters:
570eed5f15bSPeter Brune +  snes - the preconditioner context
571eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
572eed5f15bSPeter Brune 
573eed5f15bSPeter Brune    Options Database Key:
574eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
575eed5f15bSPeter Brune 
576eed5f15bSPeter Brune    Level: Developer
577eed5f15bSPeter Brune 
578eed5f15bSPeter Brune @*/
579eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
580eed5f15bSPeter Brune {
581eed5f15bSPeter Brune   PetscFunctionBegin;
582eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
583eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
584cac4c232SBarry Smith   PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));
585eed5f15bSPeter Brune   PetscFunctionReturn(0);
586eed5f15bSPeter Brune }
587eed5f15bSPeter Brune 
588eed5f15bSPeter Brune /*@C
589eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
590eed5f15bSPeter Brune 
591eed5f15bSPeter Brune    Collective on SNES
592eed5f15bSPeter Brune 
593eed5f15bSPeter Brune    Input Parameters:
594eed5f15bSPeter Brune +  snes - the preconditioner context
595eed5f15bSPeter Brune -  type - the type of the new preconditioner
596eed5f15bSPeter Brune 
597eed5f15bSPeter Brune    Level: Developer
598eed5f15bSPeter Brune 
599eed5f15bSPeter Brune @*/
600eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
601eed5f15bSPeter Brune {
602eed5f15bSPeter Brune   PetscFunctionBegin;
603eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
604cac4c232SBarry Smith   PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));
605eed5f15bSPeter Brune   PetscFunctionReturn(0);
606eed5f15bSPeter Brune }
6074ae48641SStefano Zampini 
608eed5f15bSPeter Brune /*@
609eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
610eed5f15bSPeter Brune 
611eed5f15bSPeter Brune    Not Collective
612eed5f15bSPeter Brune 
613d8d19677SJose E. Roman    Input Parameters:
614eed5f15bSPeter Brune +  snes - the preconditioner context
615eed5f15bSPeter Brune -  n - the number of the snes requested
616eed5f15bSPeter Brune 
617eed5f15bSPeter Brune    Output Parameters:
618eed5f15bSPeter Brune .  subsnes - the SNES requested
619eed5f15bSPeter Brune 
620eed5f15bSPeter Brune    Level: Developer
621eed5f15bSPeter Brune 
622db781477SPatrick Sanan .seealso: `SNESCompositeAddSNES()`
623eed5f15bSPeter Brune @*/
624eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
625eed5f15bSPeter Brune {
626eed5f15bSPeter Brune   PetscFunctionBegin;
627eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
628eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
629cac4c232SBarry Smith   PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));
630eed5f15bSPeter Brune   PetscFunctionReturn(0);
631eed5f15bSPeter Brune }
632eed5f15bSPeter Brune 
6336b2b7f7bSPatrick Farrell /*@
6346b2b7f7bSPatrick Farrell    SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES.
6356b2b7f7bSPatrick Farrell 
6366b2b7f7bSPatrick Farrell    Logically Collective on SNES
6376b2b7f7bSPatrick Farrell 
6386b2b7f7bSPatrick Farrell    Input Parameter:
6396b2b7f7bSPatrick Farrell    snes - the preconditioner context
6406b2b7f7bSPatrick Farrell 
6416b2b7f7bSPatrick Farrell    Output Parameter:
6426b2b7f7bSPatrick Farrell    n - the number of subsolvers
6436b2b7f7bSPatrick Farrell 
6446b2b7f7bSPatrick Farrell    Level: Developer
6456b2b7f7bSPatrick Farrell 
6466b2b7f7bSPatrick Farrell @*/
6476b2b7f7bSPatrick Farrell PetscErrorCode  SNESCompositeGetNumber(SNES snes,PetscInt *n)
6486b2b7f7bSPatrick Farrell {
6496b2b7f7bSPatrick Farrell   SNES_Composite     *jac;
6506b2b7f7bSPatrick Farrell   SNES_CompositeLink next;
6516b2b7f7bSPatrick Farrell 
6526b2b7f7bSPatrick Farrell   PetscFunctionBegin;
6536b2b7f7bSPatrick Farrell   jac  = (SNES_Composite*)snes->data;
6546b2b7f7bSPatrick Farrell   next = jac->head;
6556b2b7f7bSPatrick Farrell 
6566b2b7f7bSPatrick Farrell   *n = 0;
6576b2b7f7bSPatrick Farrell   while (next) {
6586b2b7f7bSPatrick Farrell     *n = *n + 1;
6596b2b7f7bSPatrick Farrell     next = next->next;
6606b2b7f7bSPatrick Farrell   }
6616b2b7f7bSPatrick Farrell   PetscFunctionReturn(0);
6626b2b7f7bSPatrick Farrell }
6636b2b7f7bSPatrick Farrell 
6648f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
6658f626970SPeter Brune {
6668f626970SPeter Brune   SNES_Composite     *jac;
6678f626970SPeter Brune   SNES_CompositeLink next;
6688f626970SPeter Brune   PetscInt           i;
6698f626970SPeter Brune 
6708f626970SPeter Brune   PetscFunctionBegin;
6718f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
6728f626970SPeter Brune   next = jac->head;
6738f626970SPeter Brune   for (i=0; i<n; i++) {
67428b400f6SJacob Faibussowitsch     PetscCheck(next->next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
6758f626970SPeter Brune     next = next->next;
6768f626970SPeter Brune   }
6778f626970SPeter Brune   next->dmp = dmp;
6788f626970SPeter Brune   PetscFunctionReturn(0);
6798f626970SPeter Brune }
6808f626970SPeter Brune 
6818f626970SPeter Brune /*@
6828f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
6838f626970SPeter Brune 
6848f626970SPeter Brune    Not Collective
6858f626970SPeter Brune 
686d8d19677SJose E. Roman    Input Parameters:
6878f626970SPeter Brune +  snes - the preconditioner context
6888f626970SPeter Brune .  n - the number of the snes requested
6898f626970SPeter Brune -  dmp - the damping
6908f626970SPeter Brune 
6918f626970SPeter Brune    Level: Developer
6928f626970SPeter Brune 
693db781477SPatrick Sanan .seealso: `SNESCompositeAddSNES()`
6948f626970SPeter Brune @*/
6958f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
6968f626970SPeter Brune {
6978f626970SPeter Brune   PetscFunctionBegin;
6988f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
699cac4c232SBarry Smith   PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));
7008f626970SPeter Brune   PetscFunctionReturn(0);
7018f626970SPeter Brune }
7028f626970SPeter Brune 
70340244768SBarry Smith static PetscErrorCode SNESSolve_Composite(SNES snes)
704eed5f15bSPeter Brune {
7054ae48641SStefano Zampini   Vec              F,X,B,Y;
706eed5f15bSPeter Brune   PetscInt         i;
707b22975d2SPatrick Farrell   PetscReal        fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
70872edecb9SPeter Brune   SNESNormSchedule normtype;
709eed5f15bSPeter Brune   SNES_Composite   *comp = (SNES_Composite*)snes->data;
710eed5f15bSPeter Brune 
711eed5f15bSPeter Brune   PetscFunctionBegin;
712eed5f15bSPeter Brune   X = snes->vec_sol;
713eed5f15bSPeter Brune   F = snes->vec_func;
714eed5f15bSPeter Brune   B = snes->vec_rhs;
7154ae48641SStefano Zampini   Y = snes->vec_sol_update;
716eed5f15bSPeter Brune 
7179566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
718eed5f15bSPeter Brune   snes->iter   = 0;
719eed5f15bSPeter Brune   snes->norm   = 0.;
7200b762f1fSPatrick Farrell   comp->innerFailures = 0;
7219566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
722eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
7239566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normtype));
72472edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
725eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
7269566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes,X,F));
727eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
728eed5f15bSPeter Brune 
729a6da83ecSPatrick Farrell     if (snes->xl && snes->xu) {
7309566063dSJacob Faibussowitsch       PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
731a6da83ecSPatrick Farrell     } else {
7329566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
733a6da83ecSPatrick Farrell     }
734422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
7359566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
736eed5f15bSPeter Brune     snes->iter = 0;
737eed5f15bSPeter Brune     snes->norm = fnorm;
7389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7399566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
7409566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,0,snes->norm));
741eed5f15bSPeter Brune 
742eed5f15bSPeter Brune     /* test convergence */
743*dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes,converged ,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
744eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
745eed5f15bSPeter Brune   } else {
7469566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7479566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
7489566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,0,snes->norm));
749eed5f15bSPeter Brune   }
750eed5f15bSPeter Brune 
7514ae48641SStefano Zampini   for (i = 0; i < snes->max_its; i++) {
752eed5f15bSPeter Brune     /* Call general purpose update function */
753*dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes,update, snes->iter);
754eed5f15bSPeter Brune 
755b22975d2SPatrick Farrell     /* Copy the state before modification by application of the composite solver;
756b22975d2SPatrick Farrell        we will subtract the new state after application */
7579566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
758b22975d2SPatrick Farrell 
759eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
7609566063dSJacob Faibussowitsch       PetscCall(SNESCompositeApply_Additive(snes,X,B,F,&fnorm));
761eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
7629566063dSJacob Faibussowitsch       PetscCall(SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm));
76390a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
7649566063dSJacob Faibussowitsch       PetscCall(SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm));
7656c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
766d5a53f06SPeter Brune     if (snes->reason < 0) break;
767d5a53f06SPeter Brune 
768b22975d2SPatrick Farrell     /* Compute the solution update for convergence testing */
7699566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1.0, X));
770b22975d2SPatrick Farrell 
771eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
7729566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes,X,F));
773b22975d2SPatrick Farrell 
774a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
7759566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7769566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(Y, NORM_2, &snorm));
7779566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
7789566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7799566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(Y, NORM_2, &snorm));
780a6da83ecSPatrick Farrell       } else {
7819566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(F, NORM_2, &fnorm));
7829566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7839566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(Y, NORM_2, &snorm));
784b22975d2SPatrick Farrell 
7859566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(F, NORM_2, &fnorm));
7869566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7879566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(Y, NORM_2, &snorm));
788a6da83ecSPatrick Farrell       }
789422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
790b22975d2SPatrick Farrell     } else if (normtype == SNES_NORM_ALWAYS) {
7919566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7929566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(Y, NORM_2, &snorm));
7939566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7949566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(Y, NORM_2, &snorm));
795eed5f15bSPeter Brune     }
796eed5f15bSPeter Brune     /* Monitor convergence */
7979566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
798eed5f15bSPeter Brune     snes->iter = i+1;
799eed5f15bSPeter Brune     snes->norm = fnorm;
800c1e67a49SFande Kong     snes->xnorm = xnorm;
801c1e67a49SFande Kong     snes->ynorm = snorm;
8029566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
8039566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
8049566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
805eed5f15bSPeter Brune     /* Test for convergence */
806*dbbe0bcdSBarry Smith     if (normtype == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes,converged ,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);
807eed5f15bSPeter Brune     if (snes->reason) break;
808eed5f15bSPeter Brune   }
80972edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS) {
810eed5f15bSPeter Brune     if (i == snes->max_its) {
81163a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",snes->max_its));
812eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
813eed5f15bSPeter Brune     }
814eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
815eed5f15bSPeter Brune   PetscFunctionReturn(0);
816eed5f15bSPeter Brune }
817eed5f15bSPeter Brune 
818eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
819eed5f15bSPeter Brune 
820eed5f15bSPeter Brune /*MC
821eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
822eed5f15bSPeter Brune 
823eed5f15bSPeter Brune    Options Database Keys:
824eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
825eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
826eed5f15bSPeter Brune 
827eed5f15bSPeter Brune    Level: intermediate
828eed5f15bSPeter Brune 
829db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNESSetType()`, `SNESType`, `SNES`,
830db781477SPatrick Sanan           `SNESSHELL`, `SNESCompositeSetType()`, `SNESCompositeSpecialSetAlpha()`, `SNESCompositeAddSNES()`,
831db781477SPatrick Sanan           `SNESCompositeGetSNES()`
832eed5f15bSPeter Brune 
8334f02bc6aSBarry Smith    References:
834606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
8354f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
8364f02bc6aSBarry Smith 
837eed5f15bSPeter Brune M*/
838eed5f15bSPeter Brune 
839eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
840eed5f15bSPeter Brune {
841eed5f15bSPeter Brune   SNES_Composite *jac;
842eed5f15bSPeter Brune 
843eed5f15bSPeter Brune   PetscFunctionBegin;
8449566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&jac));
845eed5f15bSPeter Brune 
846eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
847eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
848eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
849eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
850eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
851eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
852eed5f15bSPeter Brune 
853d8d34be6SBarry Smith   snes->usesksp        = PETSC_FALSE;
854d8d34be6SBarry Smith 
8554fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8564fc747eaSLawrence Mitchell 
857eed5f15bSPeter Brune   snes->data = (void*)jac;
85890a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
85990a8ba9bSPeter Brune   jac->Fes   = NULL;
86090a8ba9bSPeter Brune   jac->Xes   = NULL;
8619c2f3473SPeter Brune   jac->fnorms = NULL;
86290a8ba9bSPeter Brune   jac->nsnes = 0;
8639e5d0892SLisandro Dalcin   jac->head  = NULL;
8645e806d2eSPeter Brune   jac->stol  = 0.1;
8655e806d2eSPeter Brune   jac->rtol  = 1.1;
866eed5f15bSPeter Brune 
86790a8ba9bSPeter Brune   jac->h     = NULL;
86890a8ba9bSPeter Brune   jac->s     = NULL;
86990a8ba9bSPeter Brune   jac->beta  = NULL;
87090a8ba9bSPeter Brune   jac->work  = NULL;
87190a8ba9bSPeter Brune   jac->rwork = NULL;
87290a8ba9bSPeter Brune 
8739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite));
8749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite));
8759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite));
8769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite));
877eed5f15bSPeter Brune   PetscFunctionReturn(0);
878eed5f15bSPeter Brune }
879