xref: /petsc/src/snes/impls/composite/snescomposite.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
1eed5f15bSPeter Brune /*
2eed5f15bSPeter Brune       Defines a SNES that can consist of a collection of SNESes
3eed5f15bSPeter Brune */
4af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
590a8ba9bSPeter Brune #include <petscblaslapack.h>
6eed5f15bSPeter Brune 
79e5d0892SLisandro Dalcin const char *const SNESCompositeTypes[] = {"ADDITIVE", "MULTIPLICATIVE", "ADDITIVEOPTIMAL", "SNESCompositeType", "SNES_COMPOSITE", NULL};
8eed5f15bSPeter Brune 
9eed5f15bSPeter Brune typedef struct _SNES_CompositeLink *SNES_CompositeLink;
10eed5f15bSPeter Brune struct _SNES_CompositeLink {
11eed5f15bSPeter Brune   SNES               snes;
128f626970SPeter Brune   PetscReal          dmp;
1390a8ba9bSPeter Brune   Vec                X;
14eed5f15bSPeter Brune   SNES_CompositeLink next;
15eed5f15bSPeter Brune   SNES_CompositeLink previous;
16eed5f15bSPeter Brune };
17eed5f15bSPeter Brune 
18eed5f15bSPeter Brune typedef struct {
19eed5f15bSPeter Brune   SNES_CompositeLink head;
2090a8ba9bSPeter Brune   PetscInt           nsnes;
21eed5f15bSPeter Brune   SNESCompositeType  type;
22eed5f15bSPeter Brune   Vec                Xorig;
230b762f1fSPatrick Farrell   PetscInt           innerFailures; /* the number of inner failures we've seen */
2490a8ba9bSPeter Brune 
2590a8ba9bSPeter Brune   /* context for ADDITIVEOPTIMAL */
2690a8ba9bSPeter Brune   Vec         *Xes, *Fes; /* solution and residual vectors for the subsolvers */
279c2f3473SPeter Brune   PetscReal   *fnorms;    /* norms of the solutions */
2890a8ba9bSPeter Brune   PetscScalar *h;         /* the matrix formed as q_ij = (rdot_i, rdot_j) */
299c2f3473SPeter Brune   PetscScalar *g;         /* the dotproducts of the previous function with the candidate functions */
309c2f3473SPeter Brune   PetscBLASInt n;         /* matrix dimension -- nsnes */
3190a8ba9bSPeter Brune   PetscBLASInt nrhs;      /* the number of right hand sides */
3290a8ba9bSPeter Brune   PetscBLASInt lda;       /* the padded matrix dimension */
3390a8ba9bSPeter Brune   PetscBLASInt ldb;       /* the padded vector dimension */
3490a8ba9bSPeter Brune   PetscReal   *s;         /* the singular values */
359c2f3473SPeter Brune   PetscScalar *beta;      /* the RHS and combination */
3690a8ba9bSPeter Brune   PetscReal    rcond;     /* the exit condition */
3790a8ba9bSPeter Brune   PetscBLASInt rank;      /* the effective rank */
3890a8ba9bSPeter Brune   PetscScalar *work;      /* the work vector */
3990a8ba9bSPeter Brune   PetscReal   *rwork;     /* the real work vector used for complex */
4090a8ba9bSPeter Brune   PetscBLASInt lwork;     /* the size of the work vector */
4190a8ba9bSPeter Brune   PetscBLASInt info;      /* the output condition */
4290a8ba9bSPeter Brune 
435e806d2eSPeter Brune   PetscReal rtol; /* restart tolerance for accepting the combination */
445e806d2eSPeter Brune   PetscReal stol; /* restart tolerance for the combination */
45eed5f15bSPeter Brune } SNES_Composite;
46eed5f15bSPeter Brune 
47d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm)
48d71ae5a4SJacob Faibussowitsch {
49eed5f15bSPeter Brune   SNES_Composite     *jac  = (SNES_Composite *)snes->data;
50eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
51eed5f15bSPeter Brune   Vec                 FSub;
52d5a53f06SPeter Brune   SNESConvergedReason reason;
53eed5f15bSPeter Brune 
54eed5f15bSPeter Brune   PetscFunctionBegin;
5528b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
5648a46eb9SPierre Jolivet   if (snes->normschedule == SNES_NORM_ALWAYS) PetscCall(SNESSetInitialFunction(next->snes, F));
579566063dSJacob Faibussowitsch   PetscCall(SNESSolve(next->snes, B, X));
589566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(next->snes, &reason));
59d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
600b762f1fSPatrick Farrell     jac->innerFailures++;
610b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
62d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
633ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
64d5a53f06SPeter Brune     }
650b762f1fSPatrick Farrell   }
66eed5f15bSPeter Brune 
67eed5f15bSPeter Brune   while (next->next) {
68eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
69efd4aadfSBarry Smith     if (next->snes->npcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
709566063dSJacob Faibussowitsch       PetscCall(SNESGetFunction(next->snes, &FSub, NULL, NULL));
71eed5f15bSPeter Brune       next = next->next;
729566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(next->snes, FSub));
73eed5f15bSPeter Brune     } else {
74eed5f15bSPeter Brune       next = next->next;
75eed5f15bSPeter Brune     }
769566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next->snes, B, X));
779566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next->snes, &reason));
78d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
790b762f1fSPatrick Farrell       jac->innerFailures++;
800b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
81d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
823ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
83d5a53f06SPeter Brune       }
84eed5f15bSPeter Brune     }
850b762f1fSPatrick Farrell   }
86efd4aadfSBarry Smith   if (next->snes->npcside == PC_RIGHT) {
879566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(next->snes, &FSub, NULL, NULL));
889566063dSJacob Faibussowitsch     PetscCall(VecCopy(FSub, F));
89d5a53f06SPeter Brune     if (fnorm) {
9063cdc2bdSPatrick Farrell       if (snes->xl && snes->xu) {
919566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
9263cdc2bdSPatrick Farrell       } else {
939566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, fnorm));
9463cdc2bdSPatrick Farrell       }
95422a814eSBarry Smith       SNESCheckFunctionNorm(snes, *fnorm);
96d5a53f06SPeter Brune     }
9772edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
989566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F));
99d5a53f06SPeter Brune     if (fnorm) {
100a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
1019566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
102a6da83ecSPatrick Farrell       } else {
1039566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, fnorm));
104a6da83ecSPatrick Farrell       }
105422a814eSBarry Smith       SNESCheckFunctionNorm(snes, *fnorm);
106d5a53f06SPeter Brune     }
107eed5f15bSPeter Brune   }
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109eed5f15bSPeter Brune }
110eed5f15bSPeter Brune 
111d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeApply_Additive(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm)
112d71ae5a4SJacob Faibussowitsch {
113eed5f15bSPeter Brune   SNES_Composite     *jac  = (SNES_Composite *)snes->data;
114eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
115eed5f15bSPeter Brune   Vec                 Y, Xorig;
116d5a53f06SPeter Brune   SNESConvergedReason reason;
117eed5f15bSPeter Brune 
118eed5f15bSPeter Brune   PetscFunctionBegin;
119eed5f15bSPeter Brune   Y = snes->vec_sol_update;
1209566063dSJacob Faibussowitsch   if (!jac->Xorig) PetscCall(VecDuplicate(X, &jac->Xorig));
121eed5f15bSPeter Brune   Xorig = jac->Xorig;
1229566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Xorig));
12328b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
12472edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
1259566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next->snes, F));
126eed5f15bSPeter Brune     while (next->next) {
127eed5f15bSPeter Brune       next = next->next;
1289566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(next->snes, F));
129eed5f15bSPeter Brune     }
130eed5f15bSPeter Brune   }
131eed5f15bSPeter Brune   next = jac->head;
1329566063dSJacob Faibussowitsch   PetscCall(VecCopy(Xorig, Y));
1339566063dSJacob Faibussowitsch   PetscCall(SNESSolve(next->snes, B, Y));
1349566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(next->snes, &reason));
135d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1360b762f1fSPatrick Farrell     jac->innerFailures++;
1370b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
138d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1393ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
140d5a53f06SPeter Brune     }
1410b762f1fSPatrick Farrell   }
1429566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y, -1.0, Xorig));
1439566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, next->dmp, Y));
1448f626970SPeter Brune   while (next->next) {
1458f626970SPeter Brune     next = next->next;
1469566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xorig, Y));
1479566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next->snes, B, Y));
1489566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next->snes, &reason));
149d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1500b762f1fSPatrick Farrell       jac->innerFailures++;
1510b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
152d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1533ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
154d5a53f06SPeter Brune       }
1550b762f1fSPatrick Farrell     }
1569566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, -1.0, Xorig));
1579566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X, next->dmp, Y));
158eed5f15bSPeter Brune   }
15972edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
1609566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F));
161a6da83ecSPatrick Farrell     if (fnorm) {
162a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
1639566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
164a6da83ecSPatrick Farrell       } else {
1659566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, fnorm));
166a6da83ecSPatrick Farrell       }
167422a814eSBarry Smith       SNESCheckFunctionNorm(snes, *fnorm);
168a6da83ecSPatrick Farrell     }
169eed5f15bSPeter Brune   }
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171eed5f15bSPeter Brune }
17290a8ba9bSPeter Brune 
17390a8ba9bSPeter Brune /* approximately solve the overdetermined system:
17490a8ba9bSPeter Brune 
17590a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
17690a8ba9bSPeter Brune  \alpha_i                      = 1
17790a8ba9bSPeter Brune 
17890a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
17990a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
18090a8ba9bSPeter Brune 
18190a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
18290a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
18390a8ba9bSPeter Brune  */
184d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm)
185d71ae5a4SJacob Faibussowitsch {
18690a8ba9bSPeter Brune   SNES_Composite     *jac  = (SNES_Composite *)snes->data;
18790a8ba9bSPeter Brune   SNES_CompositeLink  next = jac->head;
18890a8ba9bSPeter Brune   Vec                *Xes = jac->Xes, *Fes = jac->Fes;
18990a8ba9bSPeter Brune   PetscInt            i, j;
1905e806d2eSPeter Brune   PetscScalar         tot, total, ftf;
1915e806d2eSPeter Brune   PetscReal           min_fnorm;
1925e806d2eSPeter Brune   PetscInt            min_i;
193d5a53f06SPeter Brune   SNESConvergedReason reason;
19490a8ba9bSPeter Brune 
19590a8ba9bSPeter Brune   PetscFunctionBegin;
19628b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
19758bdfa14SPeter Brune 
19858bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
19958bdfa14SPeter Brune     next = jac->head;
2009566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next->snes, F));
20158bdfa14SPeter Brune     while (next->next) {
20258bdfa14SPeter Brune       next = next->next;
2039566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(next->snes, F));
20458bdfa14SPeter Brune     }
20558bdfa14SPeter Brune   }
20658bdfa14SPeter Brune 
20790a8ba9bSPeter Brune   next = jac->head;
20890a8ba9bSPeter Brune   i    = 0;
2099566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Xes[i]));
2109566063dSJacob Faibussowitsch   PetscCall(SNESSolve(next->snes, B, Xes[i]));
2119566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(next->snes, &reason));
212d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2130b762f1fSPatrick Farrell     jac->innerFailures++;
2140b762f1fSPatrick Farrell     if (jac->innerFailures >= snes->maxFailures) {
215d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
2163ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
217d5a53f06SPeter Brune     }
2180b762f1fSPatrick Farrell   }
21990a8ba9bSPeter Brune   while (next->next) {
22090a8ba9bSPeter Brune     i++;
22190a8ba9bSPeter Brune     next = next->next;
2229566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xes[i]));
2239566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next->snes, B, Xes[i]));
2249566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next->snes, &reason));
225d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2260b762f1fSPatrick Farrell       jac->innerFailures++;
2270b762f1fSPatrick Farrell       if (jac->innerFailures >= snes->maxFailures) {
228d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2293ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
230d5a53f06SPeter Brune       }
23190a8ba9bSPeter Brune     }
2320b762f1fSPatrick Farrell   }
23390a8ba9bSPeter Brune 
23490a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
23590a8ba9bSPeter Brune   for (i = 0; i < jac->n; i++) {
23648a46eb9SPierre Jolivet     for (j = 0; j < i + 1; j++) PetscCall(VecDotBegin(Fes[i], Fes[j], &jac->h[i + j * jac->n]));
2379566063dSJacob Faibussowitsch     PetscCall(VecDotBegin(Fes[i], F, &jac->g[i]));
23890a8ba9bSPeter Brune   }
2395e806d2eSPeter Brune 
24090a8ba9bSPeter Brune   for (i = 0; i < jac->n; i++) {
24190a8ba9bSPeter Brune     for (j = 0; j < i + 1; j++) {
2429566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(Fes[i], Fes[j], &jac->h[i + j * jac->n]));
2439c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j * jac->n]));
24490a8ba9bSPeter Brune     }
2459566063dSJacob Faibussowitsch     PetscCall(VecDotEnd(Fes[i], F, &jac->g[i]));
24690a8ba9bSPeter Brune   }
24790a8ba9bSPeter Brune 
2489c2f3473SPeter Brune   ftf = (*fnorm) * (*fnorm);
24990a8ba9bSPeter Brune 
25090a8ba9bSPeter Brune   for (i = 0; i < jac->n; i++) {
251ad540459SPierre Jolivet     for (j = i + 1; j < jac->n; j++) jac->h[i + j * jac->n] = jac->h[j + i * jac->n];
25290a8ba9bSPeter Brune   }
25390a8ba9bSPeter Brune 
25490a8ba9bSPeter Brune   for (i = 0; i < jac->n; i++) {
255ad540459SPierre Jolivet     for (j = 0; j < jac->n; j++) jac->h[i + j * jac->n] = jac->h[i + j * jac->n] - jac->g[j] - jac->g[i] + ftf;
2569c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2579c2f3473SPeter Brune   }
25890a8ba9bSPeter Brune 
25990a8ba9bSPeter Brune   jac->info  = 0;
26090a8ba9bSPeter Brune   jac->rcond = -1.;
2619566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
26290a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
263792fecdfSBarry 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));
26490a8ba9bSPeter Brune #else
265792fecdfSBarry 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));
26690a8ba9bSPeter Brune #endif
2679566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
26808401ef6SPierre Jolivet   PetscCheck(jac->info >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "Bad argument to GELSS");
26908401ef6SPierre Jolivet   PetscCheck(jac->info <= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD failed to converge");
2709c2f3473SPeter Brune   tot   = 0.;
2715e806d2eSPeter Brune   total = 0.;
27290a8ba9bSPeter Brune   for (i = 0; i < jac->n; i++) {
27308401ef6SPierre Jolivet     PetscCheck(!snes->errorifnotconverged || !PetscIsInfOrNanScalar(jac->beta[i]), PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD generated inconsistent output");
27463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "%" PetscInt_FMT ": %g\n", i, (double)PetscRealPart(jac->beta[i])));
2759c2f3473SPeter Brune     tot += jac->beta[i];
2765e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
27790a8ba9bSPeter Brune   }
2789566063dSJacob Faibussowitsch   PetscCall(VecScale(X, (1. - tot)));
2799566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, jac->n, jac->beta, Xes));
2809566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, F));
281a6da83ecSPatrick Farrell 
282a6da83ecSPatrick Farrell   if (snes->xl && snes->xu) {
2839566063dSJacob Faibussowitsch     PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
284a6da83ecSPatrick Farrell   } else {
2859566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, fnorm));
286a6da83ecSPatrick Farrell   }
28790a8ba9bSPeter Brune 
2885e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
2895e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
2905e806d2eSPeter Brune   min_i     = 0;
2919c2f3473SPeter Brune   for (i = 0; i < jac->n; i++) {
2925e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
2935e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
2945e806d2eSPeter Brune       min_i     = i;
2959c2f3473SPeter Brune     }
2969c2f3473SPeter Brune   }
2975e806d2eSPeter Brune 
2985e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
2995e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm * jac->rtol < *fnorm) {
3009566063dSJacob Faibussowitsch     PetscCall(VecCopy(jac->Xes[min_i], X));
3019566063dSJacob Faibussowitsch     PetscCall(VecCopy(jac->Fes[min_i], F));
3025e806d2eSPeter Brune     *fnorm = min_fnorm;
3035e806d2eSPeter Brune   }
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30590a8ba9bSPeter Brune }
30690a8ba9bSPeter Brune 
307d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_Composite(SNES snes)
308d71ae5a4SJacob Faibussowitsch {
309dd515a93SPeter Brune   DM                 dm;
310eed5f15bSPeter Brune   SNES_Composite    *jac  = (SNES_Composite *)snes->data;
311eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
31290a8ba9bSPeter Brune   PetscInt           n    = 0, i;
31390a8ba9bSPeter Brune   Vec                F;
314eed5f15bSPeter Brune 
315eed5f15bSPeter Brune   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
31763cdc2bdSPatrick Farrell 
31863cdc2bdSPatrick Farrell   if (snes->ops->computevariablebounds) {
31963cdc2bdSPatrick Farrell     /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */
3209566063dSJacob Faibussowitsch     if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3219566063dSJacob Faibussowitsch     if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
322dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
32363cdc2bdSPatrick Farrell   }
32463cdc2bdSPatrick Farrell 
325eed5f15bSPeter Brune   while (next) {
32690a8ba9bSPeter Brune     n++;
3279566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(next->snes, dm));
3289566063dSJacob Faibussowitsch     PetscCall(SNESSetApplicationContext(next->snes, snes->user));
32963cdc2bdSPatrick Farrell     if (snes->xl && snes->xu) {
33063cdc2bdSPatrick Farrell       if (snes->ops->computevariablebounds) {
3319566063dSJacob Faibussowitsch         PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds));
33263cdc2bdSPatrick Farrell       } else {
3339566063dSJacob Faibussowitsch         PetscCall(SNESVISetVariableBounds(next->snes, snes->xl, snes->xu));
33463cdc2bdSPatrick Farrell       }
33563cdc2bdSPatrick Farrell     }
33663cdc2bdSPatrick Farrell 
337eed5f15bSPeter Brune     next = next->next;
338eed5f15bSPeter Brune   }
33990a8ba9bSPeter Brune   jac->nsnes = n;
3409566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
34190a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
3429566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(F, jac->nsnes, &jac->Xes));
3439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &jac->Fes));
3449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &jac->fnorms));
34590a8ba9bSPeter Brune     next = jac->head;
34690a8ba9bSPeter Brune     i    = 0;
34790a8ba9bSPeter Brune     while (next) {
3489566063dSJacob Faibussowitsch       PetscCall(SNESGetFunction(next->snes, &F, NULL, NULL));
34990a8ba9bSPeter Brune       jac->Fes[i] = F;
3509566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)F));
35190a8ba9bSPeter Brune       next = next->next;
35290a8ba9bSPeter Brune       i++;
35390a8ba9bSPeter Brune     }
35490a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
35590a8ba9bSPeter Brune     jac->nrhs = 1;
3569c2f3473SPeter Brune     jac->lda  = jac->nsnes;
35790a8ba9bSPeter Brune     jac->ldb  = jac->nsnes;
35890a8ba9bSPeter Brune     jac->n    = jac->nsnes;
35990a8ba9bSPeter Brune 
3609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n * jac->n, &jac->h));
3619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n, &jac->beta));
3629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n, &jac->s));
3639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->n, &jac->g));
3649c2f3473SPeter Brune     jac->lwork = 12 * jac->n;
365dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
3669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->lwork, &jac->rwork));
36790a8ba9bSPeter Brune #endif
3689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(jac->lwork, &jac->work));
36990a8ba9bSPeter Brune   }
37090a8ba9bSPeter Brune 
3713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
372eed5f15bSPeter Brune }
373eed5f15bSPeter Brune 
374d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_Composite(SNES snes)
375d71ae5a4SJacob Faibussowitsch {
376eed5f15bSPeter Brune   SNES_Composite    *jac  = (SNES_Composite *)snes->data;
377eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
378eed5f15bSPeter Brune 
379eed5f15bSPeter Brune   PetscFunctionBegin;
380eed5f15bSPeter Brune   while (next) {
3819566063dSJacob Faibussowitsch     PetscCall(SNESReset(next->snes));
382eed5f15bSPeter Brune     next = next->next;
383eed5f15bSPeter Brune   }
3849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->Xorig));
3859566063dSJacob Faibussowitsch   if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Xes));
3869566063dSJacob Faibussowitsch   if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Fes));
3879566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->fnorms));
3889566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->h));
3899566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->s));
3909566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g));
3919566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->beta));
3929566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->work));
3939566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->rwork));
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
395eed5f15bSPeter Brune }
396eed5f15bSPeter Brune 
397d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_Composite(SNES snes)
398d71ae5a4SJacob Faibussowitsch {
399eed5f15bSPeter Brune   SNES_Composite    *jac  = (SNES_Composite *)snes->data;
400eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head, next_tmp;
401eed5f15bSPeter Brune 
402eed5f15bSPeter Brune   PetscFunctionBegin;
4039566063dSJacob Faibussowitsch   PetscCall(SNESReset_Composite(snes));
404eed5f15bSPeter Brune   while (next) {
4059566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&next->snes));
406eed5f15bSPeter Brune     next_tmp = next;
407eed5f15bSPeter Brune     next     = next->next;
4089566063dSJacob Faibussowitsch     PetscCall(PetscFree(next_tmp));
409eed5f15bSPeter Brune   }
4102e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", NULL));
4112e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", NULL));
4122e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", NULL));
4132e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", NULL));
4149566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
416eed5f15bSPeter Brune }
417eed5f15bSPeter Brune 
418d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_Composite(SNES snes, PetscOptionItems *PetscOptionsObject)
419d71ae5a4SJacob Faibussowitsch {
420eed5f15bSPeter Brune   SNES_Composite    *jac  = (SNES_Composite *)snes->data;
421eed5f15bSPeter Brune   PetscInt           nmax = 8, i;
422eed5f15bSPeter Brune   SNES_CompositeLink next;
423eed5f15bSPeter Brune   char              *sneses[8];
4248f626970SPeter Brune   PetscReal          dmps[8];
425eed5f15bSPeter Brune   PetscBool          flg;
426eed5f15bSPeter Brune 
427eed5f15bSPeter Brune   PetscFunctionBegin;
428d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options");
4299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_composite_type", "Type of composition", "SNESCompositeSetType", SNESCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg));
4301baa6e33SBarry Smith   if (flg) PetscCall(SNESCompositeSetType(snes, jac->type));
4319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringArray("-snes_composite_sneses", "List of composite solvers", "SNESCompositeAddSNES", sneses, &nmax, &flg));
432eed5f15bSPeter Brune   if (flg) {
433eed5f15bSPeter Brune     for (i = 0; i < nmax; i++) {
4349566063dSJacob Faibussowitsch       PetscCall(SNESCompositeAddSNES(snes, sneses[i]));
4359566063dSJacob Faibussowitsch       PetscCall(PetscFree(sneses[i])); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
436eed5f15bSPeter Brune     }
437eed5f15bSPeter Brune   }
4389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-snes_composite_damping", "Damping of the additive composite solvers", "SNESCompositeSetDamping", dmps, &nmax, &flg));
4398f626970SPeter Brune   if (flg) {
44048a46eb9SPierre Jolivet     for (i = 0; i < nmax; i++) PetscCall(SNESCompositeSetDamping(snes, i, dmps[i]));
4418f626970SPeter Brune   }
4429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_composite_stol", "Step tolerance for restart on the additive composite solvers", "", jac->stol, &jac->stol, NULL));
4439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_composite_rtol", "Residual tolerance for the additive composite solvers", "", jac->rtol, &jac->rtol, NULL));
444d0609cedSBarry Smith   PetscOptionsHeadEnd();
445eed5f15bSPeter Brune 
446eed5f15bSPeter Brune   next = jac->head;
447eed5f15bSPeter Brune   while (next) {
4489566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(next->snes));
449eed5f15bSPeter Brune     next = next->next;
450eed5f15bSPeter Brune   }
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
452eed5f15bSPeter Brune }
453eed5f15bSPeter Brune 
454d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_Composite(SNES snes, PetscViewer viewer)
455d71ae5a4SJacob Faibussowitsch {
456eed5f15bSPeter Brune   SNES_Composite    *jac  = (SNES_Composite *)snes->data;
457eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
458eed5f15bSPeter Brune   PetscBool          iascii;
459eed5f15bSPeter Brune 
460eed5f15bSPeter Brune   PetscFunctionBegin;
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
462eed5f15bSPeter Brune   if (iascii) {
4639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type - %s\n", SNESCompositeTypes[jac->type]));
4649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  SNESes on composite preconditioner follow\n"));
4659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ---------------------------------\n"));
466eed5f15bSPeter Brune   }
4671baa6e33SBarry Smith   if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer));
468eed5f15bSPeter Brune   while (next) {
4699566063dSJacob Faibussowitsch     PetscCall(SNESView(next->snes, viewer));
470eed5f15bSPeter Brune     next = next->next;
471eed5f15bSPeter Brune   }
472eed5f15bSPeter Brune   if (iascii) {
4739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
4749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ---------------------------------\n"));
475eed5f15bSPeter Brune   }
4763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
477eed5f15bSPeter Brune }
478eed5f15bSPeter Brune 
479d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetType_Composite(SNES snes, SNESCompositeType type)
480d71ae5a4SJacob Faibussowitsch {
481eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite *)snes->data;
482eed5f15bSPeter Brune 
483eed5f15bSPeter Brune   PetscFunctionBegin;
484eed5f15bSPeter Brune   jac->type = type;
4853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
486eed5f15bSPeter Brune }
487eed5f15bSPeter Brune 
488d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes, SNESType type)
489d71ae5a4SJacob Faibussowitsch {
490eed5f15bSPeter Brune   SNES_Composite    *jac;
491eed5f15bSPeter Brune   SNES_CompositeLink next, ilink;
492eed5f15bSPeter Brune   PetscInt           cnt = 0;
493eed5f15bSPeter Brune   const char        *prefix;
494d726e3a5SJed Brown   char               newprefix[20];
495eed5f15bSPeter Brune   DM                 dm;
496eed5f15bSPeter Brune 
497eed5f15bSPeter Brune   PetscFunctionBegin;
4984dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
4999e5d0892SLisandro Dalcin   ilink->next = NULL;
5009566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &ilink->snes));
5019566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
5029566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(ilink->snes, dm));
5039566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(ilink->snes, snes->abstol, snes->rtol, snes->stol, 1, snes->max_funcs));
5049566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)ilink->snes));
505eed5f15bSPeter Brune   jac  = (SNES_Composite *)snes->data;
506eed5f15bSPeter Brune   next = jac->head;
507eed5f15bSPeter Brune   if (!next) {
508eed5f15bSPeter Brune     jac->head       = ilink;
509eed5f15bSPeter Brune     ilink->previous = NULL;
510eed5f15bSPeter Brune   } else {
511eed5f15bSPeter Brune     cnt++;
512eed5f15bSPeter Brune     while (next->next) {
513eed5f15bSPeter Brune       next = next->next;
514eed5f15bSPeter Brune       cnt++;
515eed5f15bSPeter Brune     }
516eed5f15bSPeter Brune     next->next      = ilink;
517eed5f15bSPeter Brune     ilink->previous = next;
518eed5f15bSPeter Brune   }
5199566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &prefix));
5209566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(ilink->snes, prefix));
5219566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(newprefix, sizeof(newprefix), "sub_%d_", (int)cnt));
5229566063dSJacob Faibussowitsch   PetscCall(SNESAppendOptionsPrefix(ilink->snes, newprefix));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes, (PetscObject)snes, 1));
5249566063dSJacob Faibussowitsch   PetscCall(SNESSetType(ilink->snes, type));
5259566063dSJacob Faibussowitsch   PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY));
52663cdc2bdSPatrick Farrell 
5278f626970SPeter Brune   ilink->dmp = 1.0;
52890a8ba9bSPeter Brune   jac->nsnes++;
5293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
530eed5f15bSPeter Brune }
531eed5f15bSPeter Brune 
532d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes, PetscInt n, SNES *subsnes)
533d71ae5a4SJacob Faibussowitsch {
534eed5f15bSPeter Brune   SNES_Composite    *jac;
535eed5f15bSPeter Brune   SNES_CompositeLink next;
536eed5f15bSPeter Brune   PetscInt           i;
537eed5f15bSPeter Brune 
538eed5f15bSPeter Brune   PetscFunctionBegin;
539eed5f15bSPeter Brune   jac  = (SNES_Composite *)snes->data;
540eed5f15bSPeter Brune   next = jac->head;
541eed5f15bSPeter Brune   for (i = 0; i < n; i++) {
54228b400f6SJacob Faibussowitsch     PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner");
543eed5f15bSPeter Brune     next = next->next;
544eed5f15bSPeter Brune   }
545eed5f15bSPeter Brune   *subsnes = next->snes;
5463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
547eed5f15bSPeter Brune }
548eed5f15bSPeter Brune 
549e31ab4f9SPeter Brune /*@C
550eed5f15bSPeter Brune   SNESCompositeSetType - Sets the type of composite preconditioner.
551eed5f15bSPeter Brune 
552c3339decSBarry Smith   Logically Collective
553eed5f15bSPeter Brune 
554d8d19677SJose E. Roman   Input Parameters:
555eed5f15bSPeter Brune + snes - the preconditioner context
556420bcc1bSBarry Smith - type - `SNES_COMPOSITE_ADDITIVE` (default), `SNES_COMPOSITE_MULTIPLICATIVE`, or `SNES_COMPOSITE_ADDITIVEOPTIMAL`
557eed5f15bSPeter Brune 
558eed5f15bSPeter Brune   Options Database Key:
559420bcc1bSBarry Smith . -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type
560eed5f15bSPeter Brune 
561e4094ef1SJacob Faibussowitsch   Level: developer
562eed5f15bSPeter Brune 
563420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCOMPOSITE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`,
564420bcc1bSBarry Smith           `PCCompositeType`
565eed5f15bSPeter Brune @*/
566d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetType(SNES snes, SNESCompositeType type)
567d71ae5a4SJacob Faibussowitsch {
568eed5f15bSPeter Brune   PetscFunctionBegin;
569eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
570eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes, type, 2);
571cac4c232SBarry Smith   PetscTryMethod(snes, "SNESCompositeSetType_C", (SNES, SNESCompositeType), (snes, type));
5723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
573eed5f15bSPeter Brune }
574eed5f15bSPeter Brune 
575eed5f15bSPeter Brune /*@C
576f6dfbefdSBarry Smith   SNESCompositeAddSNES - Adds another `SNES` to the `SNESCOMPOSITE`
577eed5f15bSPeter Brune 
578c3339decSBarry Smith   Collective
579eed5f15bSPeter Brune 
580eed5f15bSPeter Brune   Input Parameters:
581420bcc1bSBarry Smith + snes - the `SNES` context of type `SNESCOMPOSITE`
582420bcc1bSBarry Smith - type - the `SNESType` of the new solver
583eed5f15bSPeter Brune 
584e4094ef1SJacob Faibussowitsch   Level: developer
585eed5f15bSPeter Brune 
586420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeGetSNES()`
587eed5f15bSPeter Brune @*/
588d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeAddSNES(SNES snes, SNESType type)
589d71ae5a4SJacob Faibussowitsch {
590eed5f15bSPeter Brune   PetscFunctionBegin;
591eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
592cac4c232SBarry Smith   PetscTryMethod(snes, "SNESCompositeAddSNES_C", (SNES, SNESType), (snes, type));
5933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
594eed5f15bSPeter Brune }
5954ae48641SStefano Zampini 
596eed5f15bSPeter Brune /*@
597420bcc1bSBarry Smith   SNESCompositeGetSNES - Gets one of the `SNES` objects in the `SNES` of `SNESType` `SNESCOMPOSITE`
598eed5f15bSPeter Brune 
599eed5f15bSPeter Brune   Not Collective
600eed5f15bSPeter Brune 
601d8d19677SJose E. Roman   Input Parameters:
602420bcc1bSBarry Smith + snes - the `SNES` context
603420bcc1bSBarry Smith - n    - the number of the composed `SNES` requested
604eed5f15bSPeter Brune 
6052fe279fdSBarry Smith   Output Parameter:
606f6dfbefdSBarry Smith . subsnes - the `SNES` requested
607eed5f15bSPeter Brune 
608e4094ef1SJacob Faibussowitsch   Level: developer
609eed5f15bSPeter Brune 
610420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetNumber()`
611eed5f15bSPeter Brune @*/
612d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetSNES(SNES snes, PetscInt n, SNES *subsnes)
613d71ae5a4SJacob Faibussowitsch {
614eed5f15bSPeter Brune   PetscFunctionBegin;
615eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
6164f572ea9SToby Isaac   PetscAssertPointer(subsnes, 3);
617cac4c232SBarry Smith   PetscUseMethod(snes, "SNESCompositeGetSNES_C", (SNES, PetscInt, SNES *), (snes, n, subsnes));
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
619eed5f15bSPeter Brune }
620eed5f15bSPeter Brune 
6216b2b7f7bSPatrick Farrell /*@
622f6dfbefdSBarry Smith   SNESCompositeGetNumber - Get the number of subsolvers in the `SNESCOMPOSITE`
6236b2b7f7bSPatrick Farrell 
624c3339decSBarry Smith   Logically Collective
6256b2b7f7bSPatrick Farrell 
6266b2b7f7bSPatrick Farrell   Input Parameter:
627420bcc1bSBarry Smith . snes - the `SNES` context
6286b2b7f7bSPatrick Farrell 
6296b2b7f7bSPatrick Farrell   Output Parameter:
6302fe279fdSBarry Smith . n - the number of subsolvers
6316b2b7f7bSPatrick Farrell 
632e4094ef1SJacob Faibussowitsch   Level: developer
6336b2b7f7bSPatrick Farrell 
634420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`
6356b2b7f7bSPatrick Farrell @*/
636d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetNumber(SNES snes, PetscInt *n)
637d71ae5a4SJacob Faibussowitsch {
6386b2b7f7bSPatrick Farrell   SNES_Composite    *jac;
6396b2b7f7bSPatrick Farrell   SNES_CompositeLink next;
6406b2b7f7bSPatrick Farrell 
6416b2b7f7bSPatrick Farrell   PetscFunctionBegin;
6426b2b7f7bSPatrick Farrell   jac  = (SNES_Composite *)snes->data;
6436b2b7f7bSPatrick Farrell   next = jac->head;
6446b2b7f7bSPatrick Farrell 
6456b2b7f7bSPatrick Farrell   *n = 0;
6466b2b7f7bSPatrick Farrell   while (next) {
6476b2b7f7bSPatrick Farrell     *n   = *n + 1;
6486b2b7f7bSPatrick Farrell     next = next->next;
6496b2b7f7bSPatrick Farrell   }
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6516b2b7f7bSPatrick Farrell }
6526b2b7f7bSPatrick Farrell 
653d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes, PetscInt n, PetscReal dmp)
654d71ae5a4SJacob Faibussowitsch {
6558f626970SPeter Brune   SNES_Composite    *jac;
6568f626970SPeter Brune   SNES_CompositeLink next;
6578f626970SPeter Brune   PetscInt           i;
6588f626970SPeter Brune 
6598f626970SPeter Brune   PetscFunctionBegin;
6608f626970SPeter Brune   jac  = (SNES_Composite *)snes->data;
6618f626970SPeter Brune   next = jac->head;
6628f626970SPeter Brune   for (i = 0; i < n; i++) {
66328b400f6SJacob Faibussowitsch     PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner");
6648f626970SPeter Brune     next = next->next;
6658f626970SPeter Brune   }
6668f626970SPeter Brune   next->dmp = dmp;
6673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6688f626970SPeter Brune }
6698f626970SPeter Brune 
6708f626970SPeter Brune /*@
671420bcc1bSBarry Smith   SNESCompositeSetDamping - Sets the damping of a subsolver when using `SNES_COMPOSITE_ADDITIVE` with a `SNES` of `SNESType` `SNESCOMPOSITE`
6728f626970SPeter Brune 
6738f626970SPeter Brune   Not Collective
6748f626970SPeter Brune 
675d8d19677SJose E. Roman   Input Parameters:
676420bcc1bSBarry Smith + snes - the `SNES` context
6772fe279fdSBarry Smith . n    - the number of the sub-`SNES` object requested
6788f626970SPeter Brune - dmp  - the damping
6798f626970SPeter Brune 
680f6dfbefdSBarry Smith   Level: intermediate
6818f626970SPeter Brune 
682420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`,
683f6dfbefdSBarry Smith           `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()`
6848f626970SPeter Brune @*/
685d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetDamping(SNES snes, PetscInt n, PetscReal dmp)
686d71ae5a4SJacob Faibussowitsch {
6878f626970SPeter Brune   PetscFunctionBegin;
6888f626970SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
689cac4c232SBarry Smith   PetscUseMethod(snes, "SNESCompositeSetDamping_C", (SNES, PetscInt, PetscReal), (snes, n, dmp));
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6918f626970SPeter Brune }
6928f626970SPeter Brune 
693d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_Composite(SNES snes)
694d71ae5a4SJacob Faibussowitsch {
6954ae48641SStefano Zampini   Vec              F, X, B, Y;
696eed5f15bSPeter Brune   PetscInt         i;
697b22975d2SPatrick Farrell   PetscReal        fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
69872edecb9SPeter Brune   SNESNormSchedule normtype;
699eed5f15bSPeter Brune   SNES_Composite  *comp = (SNES_Composite *)snes->data;
700eed5f15bSPeter Brune 
701eed5f15bSPeter Brune   PetscFunctionBegin;
702eed5f15bSPeter Brune   X = snes->vec_sol;
703eed5f15bSPeter Brune   F = snes->vec_func;
704eed5f15bSPeter Brune   B = snes->vec_rhs;
7054ae48641SStefano Zampini   Y = snes->vec_sol_update;
706eed5f15bSPeter Brune 
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
708eed5f15bSPeter Brune   snes->iter          = 0;
709eed5f15bSPeter Brune   snes->norm          = 0.;
7100b762f1fSPatrick Farrell   comp->innerFailures = 0;
7119566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
712eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
7139566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normtype));
71472edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
715eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
7169566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
717eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
718eed5f15bSPeter Brune 
719a6da83ecSPatrick Farrell     if (snes->xl && snes->xu) {
7209566063dSJacob Faibussowitsch       PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
721a6da83ecSPatrick Farrell     } else {
7229566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
723a6da83ecSPatrick Farrell     }
724422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
7259566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
726eed5f15bSPeter Brune     snes->iter = 0;
727eed5f15bSPeter Brune     snes->norm = fnorm;
7289566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7299566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
730eed5f15bSPeter Brune 
731eed5f15bSPeter Brune     /* test convergence */
732d76a863bSStefano Zampini     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
7332d157150SStefano Zampini     PetscCall(SNESMonitor(snes, 0, snes->norm));
7343ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
735eed5f15bSPeter Brune   } else {
7369566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7379566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7389566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 0, snes->norm));
739eed5f15bSPeter Brune   }
740eed5f15bSPeter Brune 
7414ae48641SStefano Zampini   for (i = 0; i < snes->max_its; i++) {
742eed5f15bSPeter Brune     /* Call general purpose update function */
743dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
744eed5f15bSPeter Brune 
745b22975d2SPatrick Farrell     /* Copy the state before modification by application of the composite solver;
746b22975d2SPatrick Farrell        we will subtract the new state after application */
7479566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
748b22975d2SPatrick Farrell 
749eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
7509566063dSJacob Faibussowitsch       PetscCall(SNESCompositeApply_Additive(snes, X, B, F, &fnorm));
751eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
7529566063dSJacob Faibussowitsch       PetscCall(SNESCompositeApply_Multiplicative(snes, X, B, F, &fnorm));
75390a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
7549566063dSJacob Faibussowitsch       PetscCall(SNESCompositeApply_AdditiveOptimal(snes, X, B, F, &fnorm));
7556c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported SNESComposite type");
756d5a53f06SPeter Brune     if (snes->reason < 0) break;
757d5a53f06SPeter Brune 
758b22975d2SPatrick Farrell     /* Compute the solution update for convergence testing */
7599566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1.0, X));
760b22975d2SPatrick Farrell 
761eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
7629566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
763b22975d2SPatrick Farrell 
764a6da83ecSPatrick Farrell       if (snes->xl && snes->xu) {
7659566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7669566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(Y, NORM_2, &snorm));
7679566063dSJacob Faibussowitsch         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
7689566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7699566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(Y, NORM_2, &snorm));
770a6da83ecSPatrick Farrell       } else {
7719566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(F, NORM_2, &fnorm));
7729566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7739566063dSJacob Faibussowitsch         PetscCall(VecNormBegin(Y, NORM_2, &snorm));
774b22975d2SPatrick Farrell 
7759566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(F, NORM_2, &fnorm));
7769566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7779566063dSJacob Faibussowitsch         PetscCall(VecNormEnd(Y, NORM_2, &snorm));
778a6da83ecSPatrick Farrell       }
779422a814eSBarry Smith       SNESCheckFunctionNorm(snes, fnorm);
780b22975d2SPatrick Farrell     } else if (normtype == SNES_NORM_ALWAYS) {
7819566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7829566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(Y, NORM_2, &snorm));
7839566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7849566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(Y, NORM_2, &snorm));
785eed5f15bSPeter Brune     }
786eed5f15bSPeter Brune     /* Monitor convergence */
7879566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
788eed5f15bSPeter Brune     snes->iter  = i + 1;
789eed5f15bSPeter Brune     snes->norm  = fnorm;
790c1e67a49SFande Kong     snes->xnorm = xnorm;
791c1e67a49SFande Kong     snes->ynorm = snorm;
7929566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7939566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
794eed5f15bSPeter Brune     /* Test for convergence */
7952d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, snorm, fnorm));
7962d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
797eed5f15bSPeter Brune     if (snes->reason) break;
798eed5f15bSPeter Brune   }
7993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
800eed5f15bSPeter Brune }
801eed5f15bSPeter Brune 
802eed5f15bSPeter Brune /*MC
803*1d27aa22SBarry Smith      SNESCOMPOSITE - Builds a nonlinear solver/preconditioner by composing together several `SNES` nonlinear solvers {cite}`bruneknepleysmithtu15`
804eed5f15bSPeter Brune 
805eed5f15bSPeter Brune    Options Database Keys:
806420bcc1bSBarry Smith +  -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type
8072fe279fdSBarry Smith -  -snes_composite_sneses - <snes0,snes1,...> list of `SNES` to compose
808eed5f15bSPeter Brune 
809eed5f15bSPeter Brune    Level: intermediate
810eed5f15bSPeter Brune 
811420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`,
812420bcc1bSBarry Smith           `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`,
813420bcc1bSBarry Smith           `SNESCompositeSetType()`, `SNESCompositeSetDamping()`, `SNESCompositeGetNumber()`, `PCCompositeType`
814eed5f15bSPeter Brune M*/
815eed5f15bSPeter Brune 
816d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
817d71ae5a4SJacob Faibussowitsch {
818eed5f15bSPeter Brune   SNES_Composite *jac;
819eed5f15bSPeter Brune 
820eed5f15bSPeter Brune   PetscFunctionBegin;
8214dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
822eed5f15bSPeter Brune 
823eed5f15bSPeter Brune   snes->ops->solve          = SNESSolve_Composite;
824eed5f15bSPeter Brune   snes->ops->setup          = SNESSetUp_Composite;
825eed5f15bSPeter Brune   snes->ops->reset          = SNESReset_Composite;
826eed5f15bSPeter Brune   snes->ops->destroy        = SNESDestroy_Composite;
827eed5f15bSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_Composite;
828eed5f15bSPeter Brune   snes->ops->view           = SNESView_Composite;
829eed5f15bSPeter Brune 
830d8d34be6SBarry Smith   snes->usesksp = PETSC_FALSE;
831d8d34be6SBarry Smith 
8324fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8334fc747eaSLawrence Mitchell 
834eed5f15bSPeter Brune   snes->data  = (void *)jac;
83590a8ba9bSPeter Brune   jac->type   = SNES_COMPOSITE_ADDITIVEOPTIMAL;
83690a8ba9bSPeter Brune   jac->Fes    = NULL;
83790a8ba9bSPeter Brune   jac->Xes    = NULL;
8389c2f3473SPeter Brune   jac->fnorms = NULL;
83990a8ba9bSPeter Brune   jac->nsnes  = 0;
8409e5d0892SLisandro Dalcin   jac->head   = NULL;
8415e806d2eSPeter Brune   jac->stol   = 0.1;
8425e806d2eSPeter Brune   jac->rtol   = 1.1;
843eed5f15bSPeter Brune 
84490a8ba9bSPeter Brune   jac->h     = NULL;
84590a8ba9bSPeter Brune   jac->s     = NULL;
84690a8ba9bSPeter Brune   jac->beta  = NULL;
84790a8ba9bSPeter Brune   jac->work  = NULL;
84890a8ba9bSPeter Brune   jac->rwork = NULL;
84990a8ba9bSPeter Brune 
8509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", SNESCompositeSetType_Composite));
8519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", SNESCompositeAddSNES_Composite));
8529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", SNESCompositeGetSNES_Composite));
8539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", SNESCompositeSetDamping_Composite));
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
855eed5f15bSPeter Brune }
856