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 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 } 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 } 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)); 32836e53afaSBarry Smith PetscCall(SNESSetJacobian(next->snes, snes->jacobian, snes->jacobian_pre, NULL, NULL)); 329*49abdd8aSBarry Smith PetscCall(SNESSetApplicationContext(next->snes, snes->ctx)); 33063cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 33163cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 3329566063dSJacob Faibussowitsch PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds)); 33363cdc2bdSPatrick Farrell } else { 3349566063dSJacob Faibussowitsch PetscCall(SNESVISetVariableBounds(next->snes, snes->xl, snes->xu)); 33563cdc2bdSPatrick Farrell } 33663cdc2bdSPatrick Farrell } 33763cdc2bdSPatrick Farrell 338eed5f15bSPeter Brune next = next->next; 339eed5f15bSPeter Brune } 34090a8ba9bSPeter Brune jac->nsnes = n; 3419566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 34290a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 3439566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(F, jac->nsnes, &jac->Xes)); 3449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &jac->Fes)); 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &jac->fnorms)); 34690a8ba9bSPeter Brune next = jac->head; 34790a8ba9bSPeter Brune i = 0; 34890a8ba9bSPeter Brune while (next) { 3499566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(next->snes, &F, NULL, NULL)); 35090a8ba9bSPeter Brune jac->Fes[i] = F; 3519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)F)); 35290a8ba9bSPeter Brune next = next->next; 35390a8ba9bSPeter Brune i++; 35490a8ba9bSPeter Brune } 35590a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 35690a8ba9bSPeter Brune jac->nrhs = 1; 3576497c311SBarry Smith jac->lda = (PetscBLASInt)jac->nsnes; 3586497c311SBarry Smith jac->ldb = (PetscBLASInt)jac->nsnes; 3596497c311SBarry Smith jac->n = (PetscBLASInt)jac->nsnes; 36090a8ba9bSPeter Brune 3616497c311SBarry Smith PetscCall(PetscMalloc4(jac->n * jac->n, &jac->h, jac->n, &jac->beta, jac->n, &jac->s, jac->n, &jac->g)); 3629c2f3473SPeter Brune jac->lwork = 12 * jac->n; 363dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 3649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->rwork)); 36590a8ba9bSPeter Brune #endif 3669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->work)); 36790a8ba9bSPeter Brune } 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 369eed5f15bSPeter Brune } 370eed5f15bSPeter Brune 371d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_Composite(SNES snes) 372d71ae5a4SJacob Faibussowitsch { 373eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 374eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 375eed5f15bSPeter Brune 376eed5f15bSPeter Brune PetscFunctionBegin; 377eed5f15bSPeter Brune while (next) { 3789566063dSJacob Faibussowitsch PetscCall(SNESReset(next->snes)); 379eed5f15bSPeter Brune next = next->next; 380eed5f15bSPeter Brune } 3819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->Xorig)); 3829566063dSJacob Faibussowitsch if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Xes)); 3839566063dSJacob Faibussowitsch if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Fes)); 3849566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->fnorms)); 3856497c311SBarry Smith PetscCall(PetscFree4(jac->h, jac->beta, jac->s, jac->g)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->work)); 3879566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->rwork)); 3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 389eed5f15bSPeter Brune } 390eed5f15bSPeter Brune 391d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_Composite(SNES snes) 392d71ae5a4SJacob Faibussowitsch { 393eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 394eed5f15bSPeter Brune SNES_CompositeLink next = jac->head, next_tmp; 395eed5f15bSPeter Brune 396eed5f15bSPeter Brune PetscFunctionBegin; 3979566063dSJacob Faibussowitsch PetscCall(SNESReset_Composite(snes)); 398eed5f15bSPeter Brune while (next) { 3999566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&next->snes)); 400eed5f15bSPeter Brune next_tmp = next; 401eed5f15bSPeter Brune next = next->next; 4029566063dSJacob Faibussowitsch PetscCall(PetscFree(next_tmp)); 403eed5f15bSPeter Brune } 4042e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", NULL)); 4052e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", NULL)); 4062e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", NULL)); 4072e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", NULL)); 4089566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 410eed5f15bSPeter Brune } 411eed5f15bSPeter Brune 412d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_Composite(SNES snes, PetscOptionItems *PetscOptionsObject) 413d71ae5a4SJacob Faibussowitsch { 414eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 415eed5f15bSPeter Brune PetscInt nmax = 8, i; 416eed5f15bSPeter Brune SNES_CompositeLink next; 417eed5f15bSPeter Brune char *sneses[8]; 4188f626970SPeter Brune PetscReal dmps[8]; 419eed5f15bSPeter Brune PetscBool flg; 420eed5f15bSPeter Brune 421eed5f15bSPeter Brune PetscFunctionBegin; 422d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options"); 4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_composite_type", "Type of composition", "SNESCompositeSetType", SNESCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg)); 4241baa6e33SBarry Smith if (flg) PetscCall(SNESCompositeSetType(snes, jac->type)); 4259566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-snes_composite_sneses", "List of composite solvers", "SNESCompositeAddSNES", sneses, &nmax, &flg)); 426eed5f15bSPeter Brune if (flg) { 427eed5f15bSPeter Brune for (i = 0; i < nmax; i++) { 4289566063dSJacob Faibussowitsch PetscCall(SNESCompositeAddSNES(snes, sneses[i])); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree(sneses[i])); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 430eed5f15bSPeter Brune } 431eed5f15bSPeter Brune } 4329566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-snes_composite_damping", "Damping of the additive composite solvers", "SNESCompositeSetDamping", dmps, &nmax, &flg)); 4338f626970SPeter Brune if (flg) { 43448a46eb9SPierre Jolivet for (i = 0; i < nmax; i++) PetscCall(SNESCompositeSetDamping(snes, i, dmps[i])); 4358f626970SPeter Brune } 4369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_composite_stol", "Step tolerance for restart on the additive composite solvers", "", jac->stol, &jac->stol, NULL)); 4379566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_composite_rtol", "Residual tolerance for the additive composite solvers", "", jac->rtol, &jac->rtol, NULL)); 438d0609cedSBarry Smith PetscOptionsHeadEnd(); 439eed5f15bSPeter Brune 440eed5f15bSPeter Brune next = jac->head; 441eed5f15bSPeter Brune while (next) { 4429566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(next->snes)); 443eed5f15bSPeter Brune next = next->next; 444eed5f15bSPeter Brune } 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 446eed5f15bSPeter Brune } 447eed5f15bSPeter Brune 448d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_Composite(SNES snes, PetscViewer viewer) 449d71ae5a4SJacob Faibussowitsch { 450eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 451eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 452eed5f15bSPeter Brune PetscBool iascii; 453eed5f15bSPeter Brune 454eed5f15bSPeter Brune PetscFunctionBegin; 4559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 456eed5f15bSPeter Brune if (iascii) { 4579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type - %s\n", SNESCompositeTypes[jac->type])); 4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SNESes on composite preconditioner follow\n")); 4599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 460eed5f15bSPeter Brune } 4611baa6e33SBarry Smith if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 462eed5f15bSPeter Brune while (next) { 4639566063dSJacob Faibussowitsch PetscCall(SNESView(next->snes, viewer)); 464eed5f15bSPeter Brune next = next->next; 465eed5f15bSPeter Brune } 466eed5f15bSPeter Brune if (iascii) { 4679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 4689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 469eed5f15bSPeter Brune } 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 471eed5f15bSPeter Brune } 472eed5f15bSPeter Brune 473d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetType_Composite(SNES snes, SNESCompositeType type) 474d71ae5a4SJacob Faibussowitsch { 475eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 476eed5f15bSPeter Brune 477eed5f15bSPeter Brune PetscFunctionBegin; 478eed5f15bSPeter Brune jac->type = type; 4793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 480eed5f15bSPeter Brune } 481eed5f15bSPeter Brune 482d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes, SNESType type) 483d71ae5a4SJacob Faibussowitsch { 484eed5f15bSPeter Brune SNES_Composite *jac; 485eed5f15bSPeter Brune SNES_CompositeLink next, ilink; 486eed5f15bSPeter Brune PetscInt cnt = 0; 487eed5f15bSPeter Brune const char *prefix; 488d726e3a5SJed Brown char newprefix[20]; 489eed5f15bSPeter Brune DM dm; 490eed5f15bSPeter Brune 491eed5f15bSPeter Brune PetscFunctionBegin; 4924dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ilink)); 4939e5d0892SLisandro Dalcin ilink->next = NULL; 4949566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &ilink->snes)); 49531f83259SBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes, (PetscObject)snes, 1)); 4969566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 4979566063dSJacob Faibussowitsch PetscCall(SNESSetDM(ilink->snes, dm)); 4989566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(ilink->snes, snes->abstol, snes->rtol, snes->stol, 1, snes->max_funcs)); 4999566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)ilink->snes)); 500eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data; 501eed5f15bSPeter Brune next = jac->head; 502eed5f15bSPeter Brune if (!next) { 503eed5f15bSPeter Brune jac->head = ilink; 504eed5f15bSPeter Brune ilink->previous = NULL; 505eed5f15bSPeter Brune } else { 506eed5f15bSPeter Brune cnt++; 507eed5f15bSPeter Brune while (next->next) { 508eed5f15bSPeter Brune next = next->next; 509eed5f15bSPeter Brune cnt++; 510eed5f15bSPeter Brune } 511eed5f15bSPeter Brune next->next = ilink; 512eed5f15bSPeter Brune ilink->previous = next; 513eed5f15bSPeter Brune } 5149566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 5159566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(ilink->snes, prefix)); 5169566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(newprefix, sizeof(newprefix), "sub_%d_", (int)cnt)); 5179566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(ilink->snes, newprefix)); 5189566063dSJacob Faibussowitsch PetscCall(SNESSetType(ilink->snes, type)); 5199566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY)); 52063cdc2bdSPatrick Farrell 5218f626970SPeter Brune ilink->dmp = 1.0; 52290a8ba9bSPeter Brune jac->nsnes++; 5233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 524eed5f15bSPeter Brune } 525eed5f15bSPeter Brune 526d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes, PetscInt n, SNES *subsnes) 527d71ae5a4SJacob Faibussowitsch { 528eed5f15bSPeter Brune SNES_Composite *jac; 529eed5f15bSPeter Brune SNES_CompositeLink next; 530eed5f15bSPeter Brune PetscInt i; 531eed5f15bSPeter Brune 532eed5f15bSPeter Brune PetscFunctionBegin; 533eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data; 534eed5f15bSPeter Brune next = jac->head; 535eed5f15bSPeter Brune for (i = 0; i < n; i++) { 53628b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 537eed5f15bSPeter Brune next = next->next; 538eed5f15bSPeter Brune } 539eed5f15bSPeter Brune *subsnes = next->snes; 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 541eed5f15bSPeter Brune } 542eed5f15bSPeter Brune 543cc4c1da9SBarry Smith /*@ 544eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 545eed5f15bSPeter Brune 546c3339decSBarry Smith Logically Collective 547eed5f15bSPeter Brune 548d8d19677SJose E. Roman Input Parameters: 549eed5f15bSPeter Brune + snes - the preconditioner context 550420bcc1bSBarry Smith - type - `SNES_COMPOSITE_ADDITIVE` (default), `SNES_COMPOSITE_MULTIPLICATIVE`, or `SNES_COMPOSITE_ADDITIVEOPTIMAL` 551eed5f15bSPeter Brune 552eed5f15bSPeter Brune Options Database Key: 553420bcc1bSBarry Smith . -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type 554eed5f15bSPeter Brune 555e4094ef1SJacob Faibussowitsch Level: developer 556eed5f15bSPeter Brune 557420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCOMPOSITE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`, 558420bcc1bSBarry Smith `PCCompositeType` 559eed5f15bSPeter Brune @*/ 560d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetType(SNES snes, SNESCompositeType type) 561d71ae5a4SJacob Faibussowitsch { 562eed5f15bSPeter Brune PetscFunctionBegin; 563eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 564eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes, type, 2); 565cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeSetType_C", (SNES, SNESCompositeType), (snes, type)); 5663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 567eed5f15bSPeter Brune } 568eed5f15bSPeter Brune 5695d83a8b1SBarry Smith /*@ 570f6dfbefdSBarry Smith SNESCompositeAddSNES - Adds another `SNES` to the `SNESCOMPOSITE` 571eed5f15bSPeter Brune 572c3339decSBarry Smith Collective 573eed5f15bSPeter Brune 574eed5f15bSPeter Brune Input Parameters: 575420bcc1bSBarry Smith + snes - the `SNES` context of type `SNESCOMPOSITE` 576420bcc1bSBarry Smith - type - the `SNESType` of the new solver 577eed5f15bSPeter Brune 578e4094ef1SJacob Faibussowitsch Level: developer 579eed5f15bSPeter Brune 580420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeGetSNES()` 581eed5f15bSPeter Brune @*/ 582d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeAddSNES(SNES snes, SNESType type) 583d71ae5a4SJacob Faibussowitsch { 584eed5f15bSPeter Brune PetscFunctionBegin; 585eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 586cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeAddSNES_C", (SNES, SNESType), (snes, type)); 5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 588eed5f15bSPeter Brune } 5894ae48641SStefano Zampini 590eed5f15bSPeter Brune /*@ 591420bcc1bSBarry Smith SNESCompositeGetSNES - Gets one of the `SNES` objects in the `SNES` of `SNESType` `SNESCOMPOSITE` 592eed5f15bSPeter Brune 593eed5f15bSPeter Brune Not Collective 594eed5f15bSPeter Brune 595d8d19677SJose E. Roman Input Parameters: 596420bcc1bSBarry Smith + snes - the `SNES` context 597420bcc1bSBarry Smith - n - the number of the composed `SNES` requested 598eed5f15bSPeter Brune 5992fe279fdSBarry Smith Output Parameter: 600f6dfbefdSBarry Smith . subsnes - the `SNES` requested 601eed5f15bSPeter Brune 602e4094ef1SJacob Faibussowitsch Level: developer 603eed5f15bSPeter Brune 604420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetNumber()` 605eed5f15bSPeter Brune @*/ 606d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetSNES(SNES snes, PetscInt n, SNES *subsnes) 607d71ae5a4SJacob Faibussowitsch { 608eed5f15bSPeter Brune PetscFunctionBegin; 609eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 6104f572ea9SToby Isaac PetscAssertPointer(subsnes, 3); 611cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeGetSNES_C", (SNES, PetscInt, SNES *), (snes, n, subsnes)); 6123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 613eed5f15bSPeter Brune } 614eed5f15bSPeter Brune 6156b2b7f7bSPatrick Farrell /*@ 616f6dfbefdSBarry Smith SNESCompositeGetNumber - Get the number of subsolvers in the `SNESCOMPOSITE` 6176b2b7f7bSPatrick Farrell 618c3339decSBarry Smith Logically Collective 6196b2b7f7bSPatrick Farrell 6206b2b7f7bSPatrick Farrell Input Parameter: 621420bcc1bSBarry Smith . snes - the `SNES` context 6226b2b7f7bSPatrick Farrell 6236b2b7f7bSPatrick Farrell Output Parameter: 6242fe279fdSBarry Smith . n - the number of subsolvers 6256b2b7f7bSPatrick Farrell 626e4094ef1SJacob Faibussowitsch Level: developer 6276b2b7f7bSPatrick Farrell 628420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()` 6296b2b7f7bSPatrick Farrell @*/ 630d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetNumber(SNES snes, PetscInt *n) 631d71ae5a4SJacob Faibussowitsch { 6326b2b7f7bSPatrick Farrell SNES_Composite *jac; 6336b2b7f7bSPatrick Farrell SNES_CompositeLink next; 6346b2b7f7bSPatrick Farrell 6356b2b7f7bSPatrick Farrell PetscFunctionBegin; 6366b2b7f7bSPatrick Farrell jac = (SNES_Composite *)snes->data; 6376b2b7f7bSPatrick Farrell next = jac->head; 6386b2b7f7bSPatrick Farrell 6396b2b7f7bSPatrick Farrell *n = 0; 6406b2b7f7bSPatrick Farrell while (next) { 6416b2b7f7bSPatrick Farrell *n = *n + 1; 6426b2b7f7bSPatrick Farrell next = next->next; 6436b2b7f7bSPatrick Farrell } 6443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6456b2b7f7bSPatrick Farrell } 6466b2b7f7bSPatrick Farrell 647d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes, PetscInt n, PetscReal dmp) 648d71ae5a4SJacob Faibussowitsch { 6498f626970SPeter Brune SNES_Composite *jac; 6508f626970SPeter Brune SNES_CompositeLink next; 6518f626970SPeter Brune PetscInt i; 6528f626970SPeter Brune 6538f626970SPeter Brune PetscFunctionBegin; 6548f626970SPeter Brune jac = (SNES_Composite *)snes->data; 6558f626970SPeter Brune next = jac->head; 6568f626970SPeter Brune for (i = 0; i < n; i++) { 65728b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 6588f626970SPeter Brune next = next->next; 6598f626970SPeter Brune } 6608f626970SPeter Brune next->dmp = dmp; 6613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6628f626970SPeter Brune } 6638f626970SPeter Brune 6648f626970SPeter Brune /*@ 665420bcc1bSBarry Smith SNESCompositeSetDamping - Sets the damping of a subsolver when using `SNES_COMPOSITE_ADDITIVE` with a `SNES` of `SNESType` `SNESCOMPOSITE` 6668f626970SPeter Brune 6678f626970SPeter Brune Not Collective 6688f626970SPeter Brune 669d8d19677SJose E. Roman Input Parameters: 670420bcc1bSBarry Smith + snes - the `SNES` context 6712fe279fdSBarry Smith . n - the number of the sub-`SNES` object requested 6728f626970SPeter Brune - dmp - the damping 6738f626970SPeter Brune 674f6dfbefdSBarry Smith Level: intermediate 6758f626970SPeter Brune 676420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 677f6dfbefdSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()` 6788f626970SPeter Brune @*/ 679d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetDamping(SNES snes, PetscInt n, PetscReal dmp) 680d71ae5a4SJacob Faibussowitsch { 6818f626970SPeter Brune PetscFunctionBegin; 6828f626970SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 683cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeSetDamping_C", (SNES, PetscInt, PetscReal), (snes, n, dmp)); 6843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6858f626970SPeter Brune } 6868f626970SPeter Brune 687d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_Composite(SNES snes) 688d71ae5a4SJacob Faibussowitsch { 6894ae48641SStefano Zampini Vec F, X, B, Y; 690eed5f15bSPeter Brune PetscInt i; 691b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 69272edecb9SPeter Brune SNESNormSchedule normtype; 693eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite *)snes->data; 694eed5f15bSPeter Brune 695eed5f15bSPeter Brune PetscFunctionBegin; 696eed5f15bSPeter Brune X = snes->vec_sol; 697eed5f15bSPeter Brune F = snes->vec_func; 698eed5f15bSPeter Brune B = snes->vec_rhs; 6994ae48641SStefano Zampini Y = snes->vec_sol_update; 700eed5f15bSPeter Brune 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 702eed5f15bSPeter Brune snes->iter = 0; 703eed5f15bSPeter Brune snes->norm = 0.; 7040b762f1fSPatrick Farrell comp->innerFailures = 0; 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 706eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 7079566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normtype)); 70872edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 709eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 7109566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 711eed5f15bSPeter Brune } 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 } 718422a814eSBarry Smith SNESCheckFunctionNorm(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 } 773422a814eSBarry Smith SNESCheckFunctionNorm(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 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