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 */
31dd8e379bSPierre Jolivet 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
SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal * fnorm)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 }
95*76c63389SBarry Smith SNESCheckFunctionDomainError(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 }
105*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, *fnorm);
106d5a53f06SPeter Brune }
107eed5f15bSPeter Brune }
1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
109eed5f15bSPeter Brune }
110eed5f15bSPeter Brune
SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal * fnorm)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 }
167*76c63389SBarry Smith SNESCheckFunctionDomainError(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 */
SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal * fnorm)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 }
27857508eceSPierre Jolivet 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 }
287*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, *fnorm);
28890a8ba9bSPeter Brune
2895e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
2905e806d2eSPeter Brune min_fnorm = jac->fnorms[0];
2915e806d2eSPeter Brune min_i = 0;
2929c2f3473SPeter Brune for (i = 0; i < jac->n; i++) {
2935e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) {
2945e806d2eSPeter Brune min_fnorm = jac->fnorms[i];
2955e806d2eSPeter Brune min_i = i;
2969c2f3473SPeter Brune }
2979c2f3473SPeter Brune }
2985e806d2eSPeter Brune
2995e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */
3005e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm * jac->rtol < *fnorm) {
3019566063dSJacob Faibussowitsch PetscCall(VecCopy(jac->Xes[min_i], X));
3029566063dSJacob Faibussowitsch PetscCall(VecCopy(jac->Fes[min_i], F));
3035e806d2eSPeter Brune *fnorm = min_fnorm;
3045e806d2eSPeter Brune }
3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
30690a8ba9bSPeter Brune }
30790a8ba9bSPeter Brune
SNESSetUp_Composite(SNES snes)308d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_Composite(SNES snes)
309d71ae5a4SJacob Faibussowitsch {
310dd515a93SPeter Brune DM dm;
311eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data;
312eed5f15bSPeter Brune SNES_CompositeLink next = jac->head;
31390a8ba9bSPeter Brune PetscInt n = 0, i;
31490a8ba9bSPeter Brune Vec F;
315eed5f15bSPeter Brune
316eed5f15bSPeter Brune PetscFunctionBegin;
3179566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
31863cdc2bdSPatrick Farrell
31963cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) {
32063cdc2bdSPatrick Farrell /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */
3219566063dSJacob Faibussowitsch if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3229566063dSJacob Faibussowitsch if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
323dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
32463cdc2bdSPatrick Farrell }
32563cdc2bdSPatrick Farrell
326eed5f15bSPeter Brune while (next) {
32790a8ba9bSPeter Brune n++;
3289566063dSJacob Faibussowitsch PetscCall(SNESSetDM(next->snes, dm));
32936e53afaSBarry Smith PetscCall(SNESSetJacobian(next->snes, snes->jacobian, snes->jacobian_pre, NULL, NULL));
33049abdd8aSBarry Smith PetscCall(SNESSetApplicationContext(next->snes, snes->ctx));
33163cdc2bdSPatrick Farrell if (snes->xl && snes->xu) {
33263cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) {
3339566063dSJacob Faibussowitsch PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds));
33463cdc2bdSPatrick Farrell } else {
3359566063dSJacob Faibussowitsch PetscCall(SNESVISetVariableBounds(next->snes, snes->xl, snes->xu));
33663cdc2bdSPatrick Farrell }
33763cdc2bdSPatrick Farrell }
33863cdc2bdSPatrick Farrell
339eed5f15bSPeter Brune next = next->next;
340eed5f15bSPeter Brune }
34190a8ba9bSPeter Brune jac->nsnes = n;
3429566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
34390a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
3449566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(F, jac->nsnes, &jac->Xes));
3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &jac->Fes));
3469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &jac->fnorms));
34790a8ba9bSPeter Brune next = jac->head;
34890a8ba9bSPeter Brune i = 0;
34990a8ba9bSPeter Brune while (next) {
3509566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(next->snes, &F, NULL, NULL));
35190a8ba9bSPeter Brune jac->Fes[i] = F;
3529566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)F));
35390a8ba9bSPeter Brune next = next->next;
35490a8ba9bSPeter Brune i++;
35590a8ba9bSPeter Brune }
35690a8ba9bSPeter Brune /* allocate the subspace direct solve area */
35790a8ba9bSPeter Brune jac->nrhs = 1;
358835f2295SStefano Zampini PetscCall(PetscBLASIntCast(jac->nsnes, &jac->lda));
359835f2295SStefano Zampini PetscCall(PetscBLASIntCast(jac->nsnes, &jac->ldb));
360835f2295SStefano Zampini PetscCall(PetscBLASIntCast(jac->nsnes, &jac->n));
36190a8ba9bSPeter Brune
3626497c311SBarry Smith PetscCall(PetscMalloc4(jac->n * jac->n, &jac->h, jac->n, &jac->beta, jac->n, &jac->s, jac->n, &jac->g));
3639c2f3473SPeter Brune jac->lwork = 12 * jac->n;
364dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
3659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->rwork));
36690a8ba9bSPeter Brune #endif
3679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->work));
36890a8ba9bSPeter Brune }
3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
370eed5f15bSPeter Brune }
371eed5f15bSPeter Brune
SNESReset_Composite(SNES snes)372d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_Composite(SNES snes)
373d71ae5a4SJacob Faibussowitsch {
374eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data;
375eed5f15bSPeter Brune SNES_CompositeLink next = jac->head;
376eed5f15bSPeter Brune
377eed5f15bSPeter Brune PetscFunctionBegin;
378eed5f15bSPeter Brune while (next) {
3799566063dSJacob Faibussowitsch PetscCall(SNESReset(next->snes));
380eed5f15bSPeter Brune next = next->next;
381eed5f15bSPeter Brune }
3829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->Xorig));
3839566063dSJacob Faibussowitsch if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Xes));
3849566063dSJacob Faibussowitsch if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Fes));
3859566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->fnorms));
3866497c311SBarry Smith PetscCall(PetscFree4(jac->h, jac->beta, jac->s, jac->g));
3879566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->work));
3889566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->rwork));
3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
390eed5f15bSPeter Brune }
391eed5f15bSPeter Brune
SNESDestroy_Composite(SNES snes)392d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_Composite(SNES snes)
393d71ae5a4SJacob Faibussowitsch {
394eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data;
395eed5f15bSPeter Brune SNES_CompositeLink next = jac->head, next_tmp;
396eed5f15bSPeter Brune
397eed5f15bSPeter Brune PetscFunctionBegin;
3989566063dSJacob Faibussowitsch PetscCall(SNESReset_Composite(snes));
399eed5f15bSPeter Brune while (next) {
4009566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&next->snes));
401eed5f15bSPeter Brune next_tmp = next;
402eed5f15bSPeter Brune next = next->next;
4039566063dSJacob Faibussowitsch PetscCall(PetscFree(next_tmp));
404eed5f15bSPeter Brune }
4052e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", NULL));
4062e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", NULL));
4072e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", NULL));
4082e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", NULL));
4099566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data));
4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
411eed5f15bSPeter Brune }
412eed5f15bSPeter Brune
SNESSetFromOptions_Composite(SNES snes,PetscOptionItems PetscOptionsObject)413ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(SNES snes, PetscOptionItems PetscOptionsObject)
414d71ae5a4SJacob Faibussowitsch {
415eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data;
416eed5f15bSPeter Brune PetscInt nmax = 8, i;
417eed5f15bSPeter Brune SNES_CompositeLink next;
418eed5f15bSPeter Brune char *sneses[8];
4198f626970SPeter Brune PetscReal dmps[8];
420eed5f15bSPeter Brune PetscBool flg;
421eed5f15bSPeter Brune
422eed5f15bSPeter Brune PetscFunctionBegin;
423d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options");
4249566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_composite_type", "Type of composition", "SNESCompositeSetType", SNESCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg));
4251baa6e33SBarry Smith if (flg) PetscCall(SNESCompositeSetType(snes, jac->type));
4269566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-snes_composite_sneses", "List of composite solvers", "SNESCompositeAddSNES", sneses, &nmax, &flg));
427eed5f15bSPeter Brune if (flg) {
428eed5f15bSPeter Brune for (i = 0; i < nmax; i++) {
4299566063dSJacob Faibussowitsch PetscCall(SNESCompositeAddSNES(snes, sneses[i]));
4309566063dSJacob Faibussowitsch PetscCall(PetscFree(sneses[i])); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
431eed5f15bSPeter Brune }
432eed5f15bSPeter Brune }
4339566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-snes_composite_damping", "Damping of the additive composite solvers", "SNESCompositeSetDamping", dmps, &nmax, &flg));
4348f626970SPeter Brune if (flg) {
43548a46eb9SPierre Jolivet for (i = 0; i < nmax; i++) PetscCall(SNESCompositeSetDamping(snes, i, dmps[i]));
4368f626970SPeter Brune }
4379566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_composite_stol", "Step tolerance for restart on the additive composite solvers", "", jac->stol, &jac->stol, NULL));
4389566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_composite_rtol", "Residual tolerance for the additive composite solvers", "", jac->rtol, &jac->rtol, NULL));
439d0609cedSBarry Smith PetscOptionsHeadEnd();
440eed5f15bSPeter Brune
441eed5f15bSPeter Brune next = jac->head;
442eed5f15bSPeter Brune while (next) {
4439566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(next->snes));
444eed5f15bSPeter Brune next = next->next;
445eed5f15bSPeter Brune }
4463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
447eed5f15bSPeter Brune }
448eed5f15bSPeter Brune
SNESView_Composite(SNES snes,PetscViewer viewer)449d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_Composite(SNES snes, PetscViewer viewer)
450d71ae5a4SJacob Faibussowitsch {
451eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data;
452eed5f15bSPeter Brune SNES_CompositeLink next = jac->head;
4539f196a02SMartin Diehl PetscBool isascii;
454eed5f15bSPeter Brune
455eed5f15bSPeter Brune PetscFunctionBegin;
4569f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
4579f196a02SMartin Diehl if (isascii) {
4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type - %s\n", SNESCompositeTypes[jac->type]));
4599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SNESes on composite preconditioner follow\n"));
4609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n"));
461eed5f15bSPeter Brune }
4629f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPushTab(viewer));
463eed5f15bSPeter Brune while (next) {
4649566063dSJacob Faibussowitsch PetscCall(SNESView(next->snes, viewer));
465eed5f15bSPeter Brune next = next->next;
466eed5f15bSPeter Brune }
4679f196a02SMartin Diehl if (isascii) {
4689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
4699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n"));
470eed5f15bSPeter Brune }
4713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
472eed5f15bSPeter Brune }
473eed5f15bSPeter Brune
SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)474d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetType_Composite(SNES snes, SNESCompositeType type)
475d71ae5a4SJacob Faibussowitsch {
476eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data;
477eed5f15bSPeter Brune
478eed5f15bSPeter Brune PetscFunctionBegin;
479eed5f15bSPeter Brune jac->type = type;
4803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
481eed5f15bSPeter Brune }
482eed5f15bSPeter Brune
SNESCompositeAddSNES_Composite(SNES snes,SNESType type)483d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes, SNESType type)
484d71ae5a4SJacob Faibussowitsch {
485eed5f15bSPeter Brune SNES_Composite *jac;
486eed5f15bSPeter Brune SNES_CompositeLink next, ilink;
487eed5f15bSPeter Brune PetscInt cnt = 0;
488eed5f15bSPeter Brune const char *prefix;
489d726e3a5SJed Brown char newprefix[20];
490eed5f15bSPeter Brune DM dm;
491eed5f15bSPeter Brune
492eed5f15bSPeter Brune PetscFunctionBegin;
4934dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ilink));
4949e5d0892SLisandro Dalcin ilink->next = NULL;
4959566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &ilink->snes));
49631f83259SBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes, (PetscObject)snes, 1));
4979566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
4989566063dSJacob Faibussowitsch PetscCall(SNESSetDM(ilink->snes, dm));
4999566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(ilink->snes, snes->abstol, snes->rtol, snes->stol, 1, snes->max_funcs));
5009566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)ilink->snes));
501eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data;
502eed5f15bSPeter Brune next = jac->head;
503eed5f15bSPeter Brune if (!next) {
504eed5f15bSPeter Brune jac->head = ilink;
505eed5f15bSPeter Brune ilink->previous = NULL;
506eed5f15bSPeter Brune } else {
507eed5f15bSPeter Brune cnt++;
508eed5f15bSPeter Brune while (next->next) {
509eed5f15bSPeter Brune next = next->next;
510eed5f15bSPeter Brune cnt++;
511eed5f15bSPeter Brune }
512eed5f15bSPeter Brune next->next = ilink;
513eed5f15bSPeter Brune ilink->previous = next;
514eed5f15bSPeter Brune }
5159566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &prefix));
5169566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(ilink->snes, prefix));
517835f2295SStefano Zampini PetscCall(PetscSNPrintf(newprefix, sizeof(newprefix), "sub_%" PetscInt_FMT "_", cnt));
5189566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(ilink->snes, newprefix));
5199566063dSJacob Faibussowitsch PetscCall(SNESSetType(ilink->snes, type));
5209566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY));
52163cdc2bdSPatrick Farrell
5228f626970SPeter Brune ilink->dmp = 1.0;
52390a8ba9bSPeter Brune jac->nsnes++;
5243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
525eed5f15bSPeter Brune }
526eed5f15bSPeter Brune
SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES * subsnes)527d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes, PetscInt n, SNES *subsnes)
528d71ae5a4SJacob Faibussowitsch {
529eed5f15bSPeter Brune SNES_Composite *jac;
530eed5f15bSPeter Brune SNES_CompositeLink next;
531eed5f15bSPeter Brune PetscInt i;
532eed5f15bSPeter Brune
533eed5f15bSPeter Brune PetscFunctionBegin;
534eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data;
535eed5f15bSPeter Brune next = jac->head;
536eed5f15bSPeter Brune for (i = 0; i < n; i++) {
53728b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner");
538eed5f15bSPeter Brune next = next->next;
539eed5f15bSPeter Brune }
540eed5f15bSPeter Brune *subsnes = next->snes;
5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
542eed5f15bSPeter Brune }
543eed5f15bSPeter Brune
544cc4c1da9SBarry Smith /*@
545eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner.
546eed5f15bSPeter Brune
547c3339decSBarry Smith Logically Collective
548eed5f15bSPeter Brune
549d8d19677SJose E. Roman Input Parameters:
550eed5f15bSPeter Brune + snes - the preconditioner context
551420bcc1bSBarry Smith - type - `SNES_COMPOSITE_ADDITIVE` (default), `SNES_COMPOSITE_MULTIPLICATIVE`, or `SNES_COMPOSITE_ADDITIVEOPTIMAL`
552eed5f15bSPeter Brune
553eed5f15bSPeter Brune Options Database Key:
554420bcc1bSBarry Smith . -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type
555eed5f15bSPeter Brune
556e4094ef1SJacob Faibussowitsch Level: developer
557eed5f15bSPeter Brune
558420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCOMPOSITE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`,
559420bcc1bSBarry Smith `PCCompositeType`
560eed5f15bSPeter Brune @*/
SNESCompositeSetType(SNES snes,SNESCompositeType type)561d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetType(SNES snes, SNESCompositeType type)
562d71ae5a4SJacob Faibussowitsch {
563eed5f15bSPeter Brune PetscFunctionBegin;
564eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
565eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes, type, 2);
566cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeSetType_C", (SNES, SNESCompositeType), (snes, type));
5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
568eed5f15bSPeter Brune }
569eed5f15bSPeter Brune
5705d83a8b1SBarry Smith /*@
571f6dfbefdSBarry Smith SNESCompositeAddSNES - Adds another `SNES` to the `SNESCOMPOSITE`
572eed5f15bSPeter Brune
573c3339decSBarry Smith Collective
574eed5f15bSPeter Brune
575eed5f15bSPeter Brune Input Parameters:
576420bcc1bSBarry Smith + snes - the `SNES` context of type `SNESCOMPOSITE`
577420bcc1bSBarry Smith - type - the `SNESType` of the new solver
578eed5f15bSPeter Brune
579e4094ef1SJacob Faibussowitsch Level: developer
580eed5f15bSPeter Brune
581420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeGetSNES()`
582eed5f15bSPeter Brune @*/
SNESCompositeAddSNES(SNES snes,SNESType type)583d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeAddSNES(SNES snes, SNESType type)
584d71ae5a4SJacob Faibussowitsch {
585eed5f15bSPeter Brune PetscFunctionBegin;
586eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
587cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeAddSNES_C", (SNES, SNESType), (snes, type));
5883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
589eed5f15bSPeter Brune }
5904ae48641SStefano Zampini
591eed5f15bSPeter Brune /*@
592420bcc1bSBarry Smith SNESCompositeGetSNES - Gets one of the `SNES` objects in the `SNES` of `SNESType` `SNESCOMPOSITE`
593eed5f15bSPeter Brune
594eed5f15bSPeter Brune Not Collective
595eed5f15bSPeter Brune
596d8d19677SJose E. Roman Input Parameters:
597420bcc1bSBarry Smith + snes - the `SNES` context
598420bcc1bSBarry Smith - n - the number of the composed `SNES` requested
599eed5f15bSPeter Brune
6002fe279fdSBarry Smith Output Parameter:
601f6dfbefdSBarry Smith . subsnes - the `SNES` requested
602eed5f15bSPeter Brune
603e4094ef1SJacob Faibussowitsch Level: developer
604eed5f15bSPeter Brune
605420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetNumber()`
606eed5f15bSPeter Brune @*/
SNESCompositeGetSNES(SNES snes,PetscInt n,SNES * subsnes)607d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetSNES(SNES snes, PetscInt n, SNES *subsnes)
608d71ae5a4SJacob Faibussowitsch {
609eed5f15bSPeter Brune PetscFunctionBegin;
610eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
6114f572ea9SToby Isaac PetscAssertPointer(subsnes, 3);
612cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeGetSNES_C", (SNES, PetscInt, SNES *), (snes, n, subsnes));
6133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
614eed5f15bSPeter Brune }
615eed5f15bSPeter Brune
6166b2b7f7bSPatrick Farrell /*@
617f6dfbefdSBarry Smith SNESCompositeGetNumber - Get the number of subsolvers in the `SNESCOMPOSITE`
6186b2b7f7bSPatrick Farrell
619c3339decSBarry Smith Logically Collective
6206b2b7f7bSPatrick Farrell
6216b2b7f7bSPatrick Farrell Input Parameter:
622420bcc1bSBarry Smith . snes - the `SNES` context
6236b2b7f7bSPatrick Farrell
6246b2b7f7bSPatrick Farrell Output Parameter:
6252fe279fdSBarry Smith . n - the number of subsolvers
6266b2b7f7bSPatrick Farrell
627e4094ef1SJacob Faibussowitsch Level: developer
6286b2b7f7bSPatrick Farrell
629420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`
6306b2b7f7bSPatrick Farrell @*/
SNESCompositeGetNumber(SNES snes,PetscInt * n)631d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetNumber(SNES snes, PetscInt *n)
632d71ae5a4SJacob Faibussowitsch {
6336b2b7f7bSPatrick Farrell SNES_Composite *jac;
6346b2b7f7bSPatrick Farrell SNES_CompositeLink next;
6356b2b7f7bSPatrick Farrell
6366b2b7f7bSPatrick Farrell PetscFunctionBegin;
6376b2b7f7bSPatrick Farrell jac = (SNES_Composite *)snes->data;
6386b2b7f7bSPatrick Farrell next = jac->head;
6396b2b7f7bSPatrick Farrell
6406b2b7f7bSPatrick Farrell *n = 0;
6416b2b7f7bSPatrick Farrell while (next) {
6426b2b7f7bSPatrick Farrell *n = *n + 1;
6436b2b7f7bSPatrick Farrell next = next->next;
6446b2b7f7bSPatrick Farrell }
6453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6466b2b7f7bSPatrick Farrell }
6476b2b7f7bSPatrick Farrell
SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)648d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes, PetscInt n, PetscReal dmp)
649d71ae5a4SJacob Faibussowitsch {
6508f626970SPeter Brune SNES_Composite *jac;
6518f626970SPeter Brune SNES_CompositeLink next;
6528f626970SPeter Brune PetscInt i;
6538f626970SPeter Brune
6548f626970SPeter Brune PetscFunctionBegin;
6558f626970SPeter Brune jac = (SNES_Composite *)snes->data;
6568f626970SPeter Brune next = jac->head;
6578f626970SPeter Brune for (i = 0; i < n; i++) {
65828b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner");
6598f626970SPeter Brune next = next->next;
6608f626970SPeter Brune }
6618f626970SPeter Brune next->dmp = dmp;
6623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6638f626970SPeter Brune }
6648f626970SPeter Brune
6658f626970SPeter Brune /*@
666420bcc1bSBarry Smith SNESCompositeSetDamping - Sets the damping of a subsolver when using `SNES_COMPOSITE_ADDITIVE` with a `SNES` of `SNESType` `SNESCOMPOSITE`
6678f626970SPeter Brune
6688f626970SPeter Brune Not Collective
6698f626970SPeter Brune
670d8d19677SJose E. Roman Input Parameters:
671420bcc1bSBarry Smith + snes - the `SNES` context
6722fe279fdSBarry Smith . n - the number of the sub-`SNES` object requested
6738f626970SPeter Brune - dmp - the damping
6748f626970SPeter Brune
675f6dfbefdSBarry Smith Level: intermediate
6768f626970SPeter Brune
677420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`,
678f6dfbefdSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()`
6798f626970SPeter Brune @*/
SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)680d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetDamping(SNES snes, PetscInt n, PetscReal dmp)
681d71ae5a4SJacob Faibussowitsch {
6828f626970SPeter Brune PetscFunctionBegin;
6838f626970SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
684cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeSetDamping_C", (SNES, PetscInt, PetscReal), (snes, n, dmp));
6853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6868f626970SPeter Brune }
6878f626970SPeter Brune
SNESSolve_Composite(SNES snes)688d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_Composite(SNES snes)
689d71ae5a4SJacob Faibussowitsch {
6904ae48641SStefano Zampini Vec F, X, B, Y;
691eed5f15bSPeter Brune PetscInt i;
692b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
69372edecb9SPeter Brune SNESNormSchedule normtype;
694eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite *)snes->data;
695eed5f15bSPeter Brune
696eed5f15bSPeter Brune PetscFunctionBegin;
697eed5f15bSPeter Brune X = snes->vec_sol;
698eed5f15bSPeter Brune F = snes->vec_func;
699eed5f15bSPeter Brune B = snes->vec_rhs;
7004ae48641SStefano Zampini Y = snes->vec_sol_update;
701eed5f15bSPeter Brune
7029566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
703eed5f15bSPeter Brune snes->iter = 0;
704eed5f15bSPeter Brune snes->norm = 0.;
7050b762f1fSPatrick Farrell comp->innerFailures = 0;
7069566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
707eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING;
7089566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normtype));
70972edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
710ac530a7eSPierre Jolivet if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
711ac530a7eSPierre Jolivet else snes->vec_func_init_set = PETSC_FALSE;
712eed5f15bSPeter Brune
713a6da83ecSPatrick Farrell if (snes->xl && snes->xu) {
7149566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
715a6da83ecSPatrick Farrell } else {
7169566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */
717a6da83ecSPatrick Farrell }
718*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
7199566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
720eed5f15bSPeter Brune snes->iter = 0;
721eed5f15bSPeter Brune snes->norm = fnorm;
7229566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7239566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
724eed5f15bSPeter Brune
725eed5f15bSPeter Brune /* test convergence */
726d76a863bSStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
7272d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, snes->norm));
7283ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
729eed5f15bSPeter Brune } else {
7309566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7319566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7329566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, snes->norm));
733eed5f15bSPeter Brune }
734eed5f15bSPeter Brune
7354ae48641SStefano Zampini for (i = 0; i < snes->max_its; i++) {
736eed5f15bSPeter Brune /* Call general purpose update function */
737dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter);
738eed5f15bSPeter Brune
739b22975d2SPatrick Farrell /* Copy the state before modification by application of the composite solver;
740b22975d2SPatrick Farrell we will subtract the new state after application */
7419566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y));
742b22975d2SPatrick Farrell
743eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) {
7449566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_Additive(snes, X, B, F, &fnorm));
745eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
7469566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_Multiplicative(snes, X, B, F, &fnorm));
74790a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
7489566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_AdditiveOptimal(snes, X, B, F, &fnorm));
7496c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported SNESComposite type");
750d5a53f06SPeter Brune if (snes->reason < 0) break;
751d5a53f06SPeter Brune
752b22975d2SPatrick Farrell /* Compute the solution update for convergence testing */
7539566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1.0, X));
754b22975d2SPatrick Farrell
755eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
7569566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F));
757b22975d2SPatrick Farrell
758a6da83ecSPatrick Farrell if (snes->xl && snes->xu) {
7599566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7609566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm));
7619566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
7629566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7639566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm));
764a6da83ecSPatrick Farrell } else {
7659566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm));
7669566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7679566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm));
768b22975d2SPatrick Farrell
7699566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm));
7709566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7719566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm));
772a6da83ecSPatrick Farrell }
773*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
774b22975d2SPatrick Farrell } else if (normtype == SNES_NORM_ALWAYS) {
7759566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm));
7769566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm));
7779566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm));
7789566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm));
779eed5f15bSPeter Brune }
780eed5f15bSPeter Brune /* Monitor convergence */
7819566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
782eed5f15bSPeter Brune snes->iter = i + 1;
783eed5f15bSPeter Brune snes->norm = fnorm;
784c1e67a49SFande Kong snes->xnorm = xnorm;
785c1e67a49SFande Kong snes->ynorm = snorm;
7869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7879566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
788eed5f15bSPeter Brune /* Test for convergence */
7892d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, snorm, fnorm));
7902d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
791eed5f15bSPeter Brune if (snes->reason) break;
792eed5f15bSPeter Brune }
7933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
794eed5f15bSPeter Brune }
795eed5f15bSPeter Brune
796eed5f15bSPeter Brune /*MC
7971d27aa22SBarry Smith SNESCOMPOSITE - Builds a nonlinear solver/preconditioner by composing together several `SNES` nonlinear solvers {cite}`bruneknepleysmithtu15`
798eed5f15bSPeter Brune
799eed5f15bSPeter Brune Options Database Keys:
800420bcc1bSBarry Smith + -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type
8012fe279fdSBarry Smith - -snes_composite_sneses - <snes0,snes1,...> list of `SNES` to compose
802eed5f15bSPeter Brune
803eed5f15bSPeter Brune Level: intermediate
804eed5f15bSPeter Brune
805420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`,
806420bcc1bSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`,
807420bcc1bSBarry Smith `SNESCompositeSetType()`, `SNESCompositeSetDamping()`, `SNESCompositeGetNumber()`, `PCCompositeType`
808eed5f15bSPeter Brune M*/
809eed5f15bSPeter Brune
SNESCreate_Composite(SNES snes)810d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
811d71ae5a4SJacob Faibussowitsch {
812eed5f15bSPeter Brune SNES_Composite *jac;
813eed5f15bSPeter Brune
814eed5f15bSPeter Brune PetscFunctionBegin;
8154dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac));
816eed5f15bSPeter Brune
817eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite;
818eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite;
819eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite;
820eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite;
821eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite;
822eed5f15bSPeter Brune snes->ops->view = SNESView_Composite;
823eed5f15bSPeter Brune
824d8d34be6SBarry Smith snes->usesksp = PETSC_FALSE;
825d8d34be6SBarry Smith
8264fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE;
8274fc747eaSLawrence Mitchell
82877e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes));
82977e5a1f9SBarry Smith
830eed5f15bSPeter Brune snes->data = (void *)jac;
83190a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL;
83290a8ba9bSPeter Brune jac->Fes = NULL;
83390a8ba9bSPeter Brune jac->Xes = NULL;
8349c2f3473SPeter Brune jac->fnorms = NULL;
83590a8ba9bSPeter Brune jac->nsnes = 0;
8369e5d0892SLisandro Dalcin jac->head = NULL;
8375e806d2eSPeter Brune jac->stol = 0.1;
8385e806d2eSPeter Brune jac->rtol = 1.1;
839eed5f15bSPeter Brune
84090a8ba9bSPeter Brune jac->h = NULL;
84190a8ba9bSPeter Brune jac->s = NULL;
84290a8ba9bSPeter Brune jac->beta = NULL;
84390a8ba9bSPeter Brune jac->work = NULL;
84490a8ba9bSPeter Brune jac->rwork = NULL;
84590a8ba9bSPeter Brune
8469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", SNESCompositeSetType_Composite));
8479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", SNESCompositeAddSNES_Composite));
8489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", SNESCompositeGetSNES_Composite));
8499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", SNESCompositeSetDamping_Composite));
8503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
851eed5f15bSPeter Brune }
852