1eed5f15bSPeter Brune 2eed5f15bSPeter Brune /* 3eed5f15bSPeter Brune Defines a SNES that can consist of a collection of SNESes 4eed5f15bSPeter Brune */ 5af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 690a8ba9bSPeter Brune #include <petscblaslapack.h> 7eed5f15bSPeter Brune 89e5d0892SLisandro Dalcin const char *const SNESCompositeTypes[] = {"ADDITIVE", "MULTIPLICATIVE", "ADDITIVEOPTIMAL", "SNESCompositeType", "SNES_COMPOSITE", NULL}; 9eed5f15bSPeter Brune 10eed5f15bSPeter Brune typedef struct _SNES_CompositeLink *SNES_CompositeLink; 11eed5f15bSPeter Brune struct _SNES_CompositeLink { 12eed5f15bSPeter Brune SNES snes; 138f626970SPeter Brune PetscReal dmp; 1490a8ba9bSPeter Brune Vec X; 15eed5f15bSPeter Brune SNES_CompositeLink next; 16eed5f15bSPeter Brune SNES_CompositeLink previous; 17eed5f15bSPeter Brune }; 18eed5f15bSPeter Brune 19eed5f15bSPeter Brune typedef struct { 20eed5f15bSPeter Brune SNES_CompositeLink head; 2190a8ba9bSPeter Brune PetscInt nsnes; 22eed5f15bSPeter Brune SNESCompositeType type; 23eed5f15bSPeter Brune Vec Xorig; 240b762f1fSPatrick Farrell PetscInt innerFailures; /* the number of inner failures we've seen */ 2590a8ba9bSPeter Brune 2690a8ba9bSPeter Brune /* context for ADDITIVEOPTIMAL */ 2790a8ba9bSPeter Brune Vec *Xes, *Fes; /* solution and residual vectors for the subsolvers */ 289c2f3473SPeter Brune PetscReal *fnorms; /* norms of the solutions */ 2990a8ba9bSPeter Brune PetscScalar *h; /* the matrix formed as q_ij = (rdot_i, rdot_j) */ 309c2f3473SPeter Brune PetscScalar *g; /* the dotproducts of the previous function with the candidate functions */ 319c2f3473SPeter Brune PetscBLASInt n; /* matrix dimension -- nsnes */ 3290a8ba9bSPeter Brune PetscBLASInt nrhs; /* the number of right hand sides */ 3390a8ba9bSPeter Brune PetscBLASInt lda; /* the padded matrix dimension */ 3490a8ba9bSPeter Brune PetscBLASInt ldb; /* the padded vector dimension */ 3590a8ba9bSPeter Brune PetscReal *s; /* the singular values */ 369c2f3473SPeter Brune PetscScalar *beta; /* the RHS and combination */ 3790a8ba9bSPeter Brune PetscReal rcond; /* the exit condition */ 3890a8ba9bSPeter Brune PetscBLASInt rank; /* the effective rank */ 3990a8ba9bSPeter Brune PetscScalar *work; /* the work vector */ 4090a8ba9bSPeter Brune PetscReal *rwork; /* the real work vector used for complex */ 4190a8ba9bSPeter Brune PetscBLASInt lwork; /* the size of the work vector */ 4290a8ba9bSPeter Brune PetscBLASInt info; /* the output condition */ 4390a8ba9bSPeter Brune 445e806d2eSPeter Brune PetscReal rtol; /* restart tolerance for accepting the combination */ 455e806d2eSPeter Brune PetscReal stol; /* restart tolerance for the combination */ 46eed5f15bSPeter Brune } SNES_Composite; 47eed5f15bSPeter Brune 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) 49d71ae5a4SJacob Faibussowitsch { 50eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 51eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 52eed5f15bSPeter Brune Vec FSub; 53d5a53f06SPeter Brune SNESConvergedReason reason; 54eed5f15bSPeter Brune 55eed5f15bSPeter Brune PetscFunctionBegin; 5628b400f6SJacob Faibussowitsch PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 5748a46eb9SPierre Jolivet if (snes->normschedule == SNES_NORM_ALWAYS) PetscCall(SNESSetInitialFunction(next->snes, F)); 589566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, X)); 599566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 60d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 610b762f1fSPatrick Farrell jac->innerFailures++; 620b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 63d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65d5a53f06SPeter Brune } 660b762f1fSPatrick Farrell } 67eed5f15bSPeter Brune 68eed5f15bSPeter Brune while (next->next) { 69eed5f15bSPeter Brune /* only copy the function over in the case where the functions correspond */ 70efd4aadfSBarry Smith if (next->snes->npcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) { 719566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(next->snes, &FSub, NULL, NULL)); 72eed5f15bSPeter Brune next = next->next; 739566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, FSub)); 74eed5f15bSPeter Brune } else { 75eed5f15bSPeter Brune next = next->next; 76eed5f15bSPeter Brune } 779566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, X)); 789566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 79d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 800b762f1fSPatrick Farrell jac->innerFailures++; 810b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 82d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84d5a53f06SPeter Brune } 85eed5f15bSPeter Brune } 860b762f1fSPatrick Farrell } 87efd4aadfSBarry Smith if (next->snes->npcside == PC_RIGHT) { 889566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(next->snes, &FSub, NULL, NULL)); 899566063dSJacob Faibussowitsch PetscCall(VecCopy(FSub, F)); 90d5a53f06SPeter Brune if (fnorm) { 9163cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 929566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 9363cdc2bdSPatrick Farrell } else { 949566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 9563cdc2bdSPatrick Farrell } 96422a814eSBarry Smith SNESCheckFunctionNorm(snes, *fnorm); 97d5a53f06SPeter Brune } 9872edecb9SPeter Brune } else if (snes->normschedule == SNES_NORM_ALWAYS) { 999566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 100d5a53f06SPeter Brune if (fnorm) { 101a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 1029566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 103a6da83ecSPatrick Farrell } else { 1049566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 105a6da83ecSPatrick Farrell } 106422a814eSBarry Smith SNESCheckFunctionNorm(snes, *fnorm); 107d5a53f06SPeter Brune } 108eed5f15bSPeter Brune } 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 110eed5f15bSPeter Brune } 111eed5f15bSPeter Brune 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeApply_Additive(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) 113d71ae5a4SJacob Faibussowitsch { 114eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 115eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 116eed5f15bSPeter Brune Vec Y, Xorig; 117d5a53f06SPeter Brune SNESConvergedReason reason; 118eed5f15bSPeter Brune 119eed5f15bSPeter Brune PetscFunctionBegin; 120eed5f15bSPeter Brune Y = snes->vec_sol_update; 1219566063dSJacob Faibussowitsch if (!jac->Xorig) PetscCall(VecDuplicate(X, &jac->Xorig)); 122eed5f15bSPeter Brune Xorig = jac->Xorig; 1239566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xorig)); 12428b400f6SJacob Faibussowitsch PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 12572edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 1269566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 127eed5f15bSPeter Brune while (next->next) { 128eed5f15bSPeter Brune next = next->next; 1299566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 130eed5f15bSPeter Brune } 131eed5f15bSPeter Brune } 132eed5f15bSPeter Brune next = jac->head; 1339566063dSJacob Faibussowitsch PetscCall(VecCopy(Xorig, Y)); 1349566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Y)); 1359566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 136d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1370b762f1fSPatrick Farrell jac->innerFailures++; 1380b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 139d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141d5a53f06SPeter Brune } 1420b762f1fSPatrick Farrell } 1439566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, -1.0, Xorig)); 1449566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, next->dmp, Y)); 1458f626970SPeter Brune while (next->next) { 1468f626970SPeter Brune next = next->next; 1479566063dSJacob Faibussowitsch PetscCall(VecCopy(Xorig, Y)); 1489566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Y)); 1499566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 150d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1510b762f1fSPatrick Farrell jac->innerFailures++; 1520b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 153d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155d5a53f06SPeter Brune } 1560b762f1fSPatrick Farrell } 1579566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, -1.0, Xorig)); 1589566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, next->dmp, Y)); 159eed5f15bSPeter Brune } 16072edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 1619566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 162a6da83ecSPatrick Farrell if (fnorm) { 163a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 1649566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 165a6da83ecSPatrick Farrell } else { 1669566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 167a6da83ecSPatrick Farrell } 168422a814eSBarry Smith SNESCheckFunctionNorm(snes, *fnorm); 169a6da83ecSPatrick Farrell } 170eed5f15bSPeter Brune } 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172eed5f15bSPeter Brune } 17390a8ba9bSPeter Brune 17490a8ba9bSPeter Brune /* approximately solve the overdetermined system: 17590a8ba9bSPeter Brune 17690a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 17790a8ba9bSPeter Brune \alpha_i = 1 17890a8ba9bSPeter Brune 17990a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 18090a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 18190a8ba9bSPeter Brune 18290a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 18390a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 18490a8ba9bSPeter Brune */ 185d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) 186d71ae5a4SJacob Faibussowitsch { 18790a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 18890a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 18990a8ba9bSPeter Brune Vec *Xes = jac->Xes, *Fes = jac->Fes; 19090a8ba9bSPeter Brune PetscInt i, j; 1915e806d2eSPeter Brune PetscScalar tot, total, ftf; 1925e806d2eSPeter Brune PetscReal min_fnorm; 1935e806d2eSPeter Brune PetscInt min_i; 194d5a53f06SPeter Brune SNESConvergedReason reason; 19590a8ba9bSPeter Brune 19690a8ba9bSPeter Brune PetscFunctionBegin; 19728b400f6SJacob Faibussowitsch PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 19858bdfa14SPeter Brune 19958bdfa14SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 20058bdfa14SPeter Brune next = jac->head; 2019566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 20258bdfa14SPeter Brune while (next->next) { 20358bdfa14SPeter Brune next = next->next; 2049566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 20558bdfa14SPeter Brune } 20658bdfa14SPeter Brune } 20758bdfa14SPeter Brune 20890a8ba9bSPeter Brune next = jac->head; 20990a8ba9bSPeter Brune i = 0; 2109566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xes[i])); 2119566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Xes[i])); 2129566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 213d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2140b762f1fSPatrick Farrell jac->innerFailures++; 2150b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 216d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218d5a53f06SPeter Brune } 2190b762f1fSPatrick Farrell } 22090a8ba9bSPeter Brune while (next->next) { 22190a8ba9bSPeter Brune i++; 22290a8ba9bSPeter Brune next = next->next; 2239566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xes[i])); 2249566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Xes[i])); 2259566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 226d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2270b762f1fSPatrick Farrell jac->innerFailures++; 2280b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 229d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 231d5a53f06SPeter Brune } 23290a8ba9bSPeter Brune } 2330b762f1fSPatrick Farrell } 23490a8ba9bSPeter Brune 23590a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 23690a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 23748a46eb9SPierre Jolivet for (j = 0; j < i + 1; j++) PetscCall(VecDotBegin(Fes[i], Fes[j], &jac->h[i + j * jac->n])); 2389566063dSJacob Faibussowitsch PetscCall(VecDotBegin(Fes[i], F, &jac->g[i])); 23990a8ba9bSPeter Brune } 2405e806d2eSPeter Brune 24190a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 24290a8ba9bSPeter Brune for (j = 0; j < i + 1; j++) { 2439566063dSJacob Faibussowitsch PetscCall(VecDotEnd(Fes[i], Fes[j], &jac->h[i + j * jac->n])); 2449c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j * jac->n])); 24590a8ba9bSPeter Brune } 2469566063dSJacob Faibussowitsch PetscCall(VecDotEnd(Fes[i], F, &jac->g[i])); 24790a8ba9bSPeter Brune } 24890a8ba9bSPeter Brune 2499c2f3473SPeter Brune ftf = (*fnorm) * (*fnorm); 25090a8ba9bSPeter Brune 25190a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 252ad540459SPierre Jolivet for (j = i + 1; j < jac->n; j++) jac->h[i + j * jac->n] = jac->h[j + i * jac->n]; 25390a8ba9bSPeter Brune } 25490a8ba9bSPeter Brune 25590a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 256ad540459SPierre 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; 2579c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 2589c2f3473SPeter Brune } 25990a8ba9bSPeter Brune 26090a8ba9bSPeter Brune jac->info = 0; 26190a8ba9bSPeter Brune jac->rcond = -1.; 2629566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 26390a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 264792fecdfSBarry 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)); 26590a8ba9bSPeter Brune #else 266792fecdfSBarry 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)); 26790a8ba9bSPeter Brune #endif 2689566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 26908401ef6SPierre Jolivet PetscCheck(jac->info >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "Bad argument to GELSS"); 27008401ef6SPierre Jolivet PetscCheck(jac->info <= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD failed to converge"); 2719c2f3473SPeter Brune tot = 0.; 2725e806d2eSPeter Brune total = 0.; 27390a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 27408401ef6SPierre Jolivet PetscCheck(!snes->errorifnotconverged || !PetscIsInfOrNanScalar(jac->beta[i]), PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD generated inconsistent output"); 27563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "%" PetscInt_FMT ": %g\n", i, (double)PetscRealPart(jac->beta[i]))); 2769c2f3473SPeter Brune tot += jac->beta[i]; 2775e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 27890a8ba9bSPeter Brune } 2799566063dSJacob Faibussowitsch PetscCall(VecScale(X, (1. - tot))); 2809566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, jac->n, jac->beta, Xes)); 2819566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 282a6da83ecSPatrick Farrell 283a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 2849566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 285a6da83ecSPatrick Farrell } else { 2869566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 287a6da83ecSPatrick Farrell } 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 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)); 3299566063dSJacob Faibussowitsch PetscCall(SNESSetApplicationContext(next->snes, snes->user)); 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; 3579c2f3473SPeter Brune jac->lda = jac->nsnes; 35890a8ba9bSPeter Brune jac->ldb = jac->nsnes; 35990a8ba9bSPeter Brune jac->n = jac->nsnes; 36090a8ba9bSPeter Brune 3619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n * jac->n, &jac->h)); 3629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n, &jac->beta)); 3639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n, &jac->s)); 3649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n, &jac->g)); 3659c2f3473SPeter Brune jac->lwork = 12 * jac->n; 366dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 3679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->rwork)); 36890a8ba9bSPeter Brune #endif 3699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->work)); 37090a8ba9bSPeter Brune } 37190a8ba9bSPeter Brune 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373eed5f15bSPeter Brune } 374eed5f15bSPeter Brune 375d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_Composite(SNES snes) 376d71ae5a4SJacob Faibussowitsch { 377eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 378eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 379eed5f15bSPeter Brune 380eed5f15bSPeter Brune PetscFunctionBegin; 381eed5f15bSPeter Brune while (next) { 3829566063dSJacob Faibussowitsch PetscCall(SNESReset(next->snes)); 383eed5f15bSPeter Brune next = next->next; 384eed5f15bSPeter Brune } 3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->Xorig)); 3869566063dSJacob Faibussowitsch if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Xes)); 3879566063dSJacob Faibussowitsch if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Fes)); 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->fnorms)); 3899566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->h)); 3909566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->s)); 3919566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->g)); 3929566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->beta)); 3939566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->work)); 3949566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->rwork)); 3953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 396eed5f15bSPeter Brune } 397eed5f15bSPeter Brune 398d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_Composite(SNES snes) 399d71ae5a4SJacob Faibussowitsch { 400eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 401eed5f15bSPeter Brune SNES_CompositeLink next = jac->head, next_tmp; 402eed5f15bSPeter Brune 403eed5f15bSPeter Brune PetscFunctionBegin; 4049566063dSJacob Faibussowitsch PetscCall(SNESReset_Composite(snes)); 405eed5f15bSPeter Brune while (next) { 4069566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&next->snes)); 407eed5f15bSPeter Brune next_tmp = next; 408eed5f15bSPeter Brune next = next->next; 4099566063dSJacob Faibussowitsch PetscCall(PetscFree(next_tmp)); 410eed5f15bSPeter Brune } 4112e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", NULL)); 4122e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", NULL)); 4132e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", NULL)); 4142e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", NULL)); 4159566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 417eed5f15bSPeter Brune } 418eed5f15bSPeter Brune 419d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_Composite(SNES snes, PetscOptionItems *PetscOptionsObject) 420d71ae5a4SJacob Faibussowitsch { 421eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 422eed5f15bSPeter Brune PetscInt nmax = 8, i; 423eed5f15bSPeter Brune SNES_CompositeLink next; 424eed5f15bSPeter Brune char *sneses[8]; 4258f626970SPeter Brune PetscReal dmps[8]; 426eed5f15bSPeter Brune PetscBool flg; 427eed5f15bSPeter Brune 428eed5f15bSPeter Brune PetscFunctionBegin; 429d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options"); 4309566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_composite_type", "Type of composition", "SNESCompositeSetType", SNESCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg)); 4311baa6e33SBarry Smith if (flg) PetscCall(SNESCompositeSetType(snes, jac->type)); 4329566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-snes_composite_sneses", "List of composite solvers", "SNESCompositeAddSNES", sneses, &nmax, &flg)); 433eed5f15bSPeter Brune if (flg) { 434eed5f15bSPeter Brune for (i = 0; i < nmax; i++) { 4359566063dSJacob Faibussowitsch PetscCall(SNESCompositeAddSNES(snes, sneses[i])); 4369566063dSJacob Faibussowitsch PetscCall(PetscFree(sneses[i])); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 437eed5f15bSPeter Brune } 438eed5f15bSPeter Brune } 4399566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-snes_composite_damping", "Damping of the additive composite solvers", "SNESCompositeSetDamping", dmps, &nmax, &flg)); 4408f626970SPeter Brune if (flg) { 44148a46eb9SPierre Jolivet for (i = 0; i < nmax; i++) PetscCall(SNESCompositeSetDamping(snes, i, dmps[i])); 4428f626970SPeter Brune } 4439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_composite_stol", "Step tolerance for restart on the additive composite solvers", "", jac->stol, &jac->stol, NULL)); 4449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_composite_rtol", "Residual tolerance for the additive composite solvers", "", jac->rtol, &jac->rtol, NULL)); 445d0609cedSBarry Smith PetscOptionsHeadEnd(); 446eed5f15bSPeter Brune 447eed5f15bSPeter Brune next = jac->head; 448eed5f15bSPeter Brune while (next) { 4499566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(next->snes)); 450eed5f15bSPeter Brune next = next->next; 451eed5f15bSPeter Brune } 4523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 453eed5f15bSPeter Brune } 454eed5f15bSPeter Brune 455d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_Composite(SNES snes, PetscViewer viewer) 456d71ae5a4SJacob Faibussowitsch { 457eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 458eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 459eed5f15bSPeter Brune PetscBool iascii; 460eed5f15bSPeter Brune 461eed5f15bSPeter Brune PetscFunctionBegin; 4629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 463eed5f15bSPeter Brune if (iascii) { 4649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type - %s\n", SNESCompositeTypes[jac->type])); 4659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SNESes on composite preconditioner follow\n")); 4669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 467eed5f15bSPeter Brune } 4681baa6e33SBarry Smith if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 469eed5f15bSPeter Brune while (next) { 4709566063dSJacob Faibussowitsch PetscCall(SNESView(next->snes, viewer)); 471eed5f15bSPeter Brune next = next->next; 472eed5f15bSPeter Brune } 473eed5f15bSPeter Brune if (iascii) { 4749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 4759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 476eed5f15bSPeter Brune } 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 478eed5f15bSPeter Brune } 479eed5f15bSPeter Brune 480d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetType_Composite(SNES snes, SNESCompositeType type) 481d71ae5a4SJacob Faibussowitsch { 482eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 483eed5f15bSPeter Brune 484eed5f15bSPeter Brune PetscFunctionBegin; 485eed5f15bSPeter Brune jac->type = type; 4863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 487eed5f15bSPeter Brune } 488eed5f15bSPeter Brune 489d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes, SNESType type) 490d71ae5a4SJacob Faibussowitsch { 491eed5f15bSPeter Brune SNES_Composite *jac; 492eed5f15bSPeter Brune SNES_CompositeLink next, ilink; 493eed5f15bSPeter Brune PetscInt cnt = 0; 494eed5f15bSPeter Brune const char *prefix; 495d726e3a5SJed Brown char newprefix[20]; 496eed5f15bSPeter Brune DM dm; 497eed5f15bSPeter Brune 498eed5f15bSPeter Brune PetscFunctionBegin; 4994dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ilink)); 5009e5d0892SLisandro Dalcin ilink->next = NULL; 5019566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &ilink->snes)); 5029566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 5039566063dSJacob Faibussowitsch PetscCall(SNESSetDM(ilink->snes, dm)); 5049566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(ilink->snes, snes->abstol, snes->rtol, snes->stol, 1, snes->max_funcs)); 5059566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)ilink->snes)); 506eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data; 507eed5f15bSPeter Brune next = jac->head; 508eed5f15bSPeter Brune if (!next) { 509eed5f15bSPeter Brune jac->head = ilink; 510eed5f15bSPeter Brune ilink->previous = NULL; 511eed5f15bSPeter Brune } else { 512eed5f15bSPeter Brune cnt++; 513eed5f15bSPeter Brune while (next->next) { 514eed5f15bSPeter Brune next = next->next; 515eed5f15bSPeter Brune cnt++; 516eed5f15bSPeter Brune } 517eed5f15bSPeter Brune next->next = ilink; 518eed5f15bSPeter Brune ilink->previous = next; 519eed5f15bSPeter Brune } 5209566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 5219566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(ilink->snes, prefix)); 5229566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(newprefix, sizeof(newprefix), "sub_%d_", (int)cnt)); 5239566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(ilink->snes, newprefix)); 5249566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes, (PetscObject)snes, 1)); 5259566063dSJacob Faibussowitsch PetscCall(SNESSetType(ilink->snes, type)); 5269566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY)); 52763cdc2bdSPatrick Farrell 5288f626970SPeter Brune ilink->dmp = 1.0; 52990a8ba9bSPeter Brune jac->nsnes++; 5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 531eed5f15bSPeter Brune } 532eed5f15bSPeter Brune 533d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes, PetscInt n, SNES *subsnes) 534d71ae5a4SJacob Faibussowitsch { 535eed5f15bSPeter Brune SNES_Composite *jac; 536eed5f15bSPeter Brune SNES_CompositeLink next; 537eed5f15bSPeter Brune PetscInt i; 538eed5f15bSPeter Brune 539eed5f15bSPeter Brune PetscFunctionBegin; 540eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data; 541eed5f15bSPeter Brune next = jac->head; 542eed5f15bSPeter Brune for (i = 0; i < n; i++) { 54328b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 544eed5f15bSPeter Brune next = next->next; 545eed5f15bSPeter Brune } 546eed5f15bSPeter Brune *subsnes = next->snes; 5473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 548eed5f15bSPeter Brune } 549eed5f15bSPeter Brune 550e31ab4f9SPeter Brune /*@C 551eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 552eed5f15bSPeter Brune 553c3339decSBarry Smith Logically Collective 554eed5f15bSPeter Brune 555d8d19677SJose E. Roman Input Parameters: 556eed5f15bSPeter Brune + snes - the preconditioner context 557f6dfbefdSBarry Smith - type - `SNES_COMPOSITE_ADDITIVE` (default), `SNES_COMPOSITE_MULTIPLICATIVE` 558eed5f15bSPeter Brune 559eed5f15bSPeter Brune Options Database Key: 560eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 561eed5f15bSPeter Brune 562*e4094ef1SJacob Faibussowitsch Level: developer 563eed5f15bSPeter Brune 564f6dfbefdSBarry Smith .seealso: `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCOMPOSITE` 565eed5f15bSPeter Brune @*/ 566d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetType(SNES snes, SNESCompositeType type) 567d71ae5a4SJacob Faibussowitsch { 568eed5f15bSPeter Brune PetscFunctionBegin; 569eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 570eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes, type, 2); 571cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeSetType_C", (SNES, SNESCompositeType), (snes, type)); 5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 573eed5f15bSPeter Brune } 574eed5f15bSPeter Brune 575eed5f15bSPeter Brune /*@C 576f6dfbefdSBarry Smith SNESCompositeAddSNES - Adds another `SNES` to the `SNESCOMPOSITE` 577eed5f15bSPeter Brune 578c3339decSBarry Smith Collective 579eed5f15bSPeter Brune 580eed5f15bSPeter Brune Input Parameters: 581f6dfbefdSBarry Smith + snes - the snes context of type `SNESCOMPOSITE` 582f6dfbefdSBarry Smith - type - the type of the new solver 583eed5f15bSPeter Brune 584*e4094ef1SJacob Faibussowitsch Level: developer 585eed5f15bSPeter Brune 586f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeGetSNES()` 587eed5f15bSPeter Brune @*/ 588d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeAddSNES(SNES snes, SNESType type) 589d71ae5a4SJacob Faibussowitsch { 590eed5f15bSPeter Brune PetscFunctionBegin; 591eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 592cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeAddSNES_C", (SNES, SNESType), (snes, type)); 5933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 594eed5f15bSPeter Brune } 5954ae48641SStefano Zampini 596eed5f15bSPeter Brune /*@ 597f6dfbefdSBarry Smith SNESCompositeGetSNES - Gets one of the `SNES` objects in the composite `SNESCOMPOSITE` 598eed5f15bSPeter Brune 599eed5f15bSPeter Brune Not Collective 600eed5f15bSPeter Brune 601d8d19677SJose E. Roman Input Parameters: 602f6dfbefdSBarry Smith + snes - the snes context 603f6dfbefdSBarry Smith - n - the number of the composed snes requested 604eed5f15bSPeter Brune 6052fe279fdSBarry Smith Output Parameter: 606f6dfbefdSBarry Smith . subsnes - the `SNES` requested 607eed5f15bSPeter Brune 608*e4094ef1SJacob Faibussowitsch Level: developer 609eed5f15bSPeter Brune 610f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()` 611eed5f15bSPeter Brune @*/ 612d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetSNES(SNES snes, PetscInt n, SNES *subsnes) 613d71ae5a4SJacob Faibussowitsch { 614eed5f15bSPeter Brune PetscFunctionBegin; 615eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 616eed5f15bSPeter Brune PetscValidPointer(subsnes, 3); 617cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeGetSNES_C", (SNES, PetscInt, SNES *), (snes, n, subsnes)); 6183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 619eed5f15bSPeter Brune } 620eed5f15bSPeter Brune 6216b2b7f7bSPatrick Farrell /*@ 622f6dfbefdSBarry Smith SNESCompositeGetNumber - Get the number of subsolvers in the `SNESCOMPOSITE` 6236b2b7f7bSPatrick Farrell 624c3339decSBarry Smith Logically Collective 6256b2b7f7bSPatrick Farrell 6266b2b7f7bSPatrick Farrell Input Parameter: 6272fe279fdSBarry Smith . snes - the snes context 6286b2b7f7bSPatrick Farrell 6296b2b7f7bSPatrick Farrell Output Parameter: 6302fe279fdSBarry Smith . n - the number of subsolvers 6316b2b7f7bSPatrick Farrell 632*e4094ef1SJacob Faibussowitsch Level: developer 6336b2b7f7bSPatrick Farrell 634f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()` 6356b2b7f7bSPatrick Farrell @*/ 636d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeGetNumber(SNES snes, PetscInt *n) 637d71ae5a4SJacob Faibussowitsch { 6386b2b7f7bSPatrick Farrell SNES_Composite *jac; 6396b2b7f7bSPatrick Farrell SNES_CompositeLink next; 6406b2b7f7bSPatrick Farrell 6416b2b7f7bSPatrick Farrell PetscFunctionBegin; 6426b2b7f7bSPatrick Farrell jac = (SNES_Composite *)snes->data; 6436b2b7f7bSPatrick Farrell next = jac->head; 6446b2b7f7bSPatrick Farrell 6456b2b7f7bSPatrick Farrell *n = 0; 6466b2b7f7bSPatrick Farrell while (next) { 6476b2b7f7bSPatrick Farrell *n = *n + 1; 6486b2b7f7bSPatrick Farrell next = next->next; 6496b2b7f7bSPatrick Farrell } 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6516b2b7f7bSPatrick Farrell } 6526b2b7f7bSPatrick Farrell 653d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes, PetscInt n, PetscReal dmp) 654d71ae5a4SJacob Faibussowitsch { 6558f626970SPeter Brune SNES_Composite *jac; 6568f626970SPeter Brune SNES_CompositeLink next; 6578f626970SPeter Brune PetscInt i; 6588f626970SPeter Brune 6598f626970SPeter Brune PetscFunctionBegin; 6608f626970SPeter Brune jac = (SNES_Composite *)snes->data; 6618f626970SPeter Brune next = jac->head; 6628f626970SPeter Brune for (i = 0; i < n; i++) { 66328b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 6648f626970SPeter Brune next = next->next; 6658f626970SPeter Brune } 6668f626970SPeter Brune next->dmp = dmp; 6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6688f626970SPeter Brune } 6698f626970SPeter Brune 6708f626970SPeter Brune /*@ 671f6dfbefdSBarry Smith SNESCompositeSetDamping - Sets the damping of a subsolver when using `SNES_COMPOSITE_ADDITIVE` `SNESCOMPOSITE` 6728f626970SPeter Brune 6738f626970SPeter Brune Not Collective 6748f626970SPeter Brune 675d8d19677SJose E. Roman Input Parameters: 676f6dfbefdSBarry Smith + snes - the snes context 6772fe279fdSBarry Smith . n - the number of the sub-`SNES` object requested 6788f626970SPeter Brune - dmp - the damping 6798f626970SPeter Brune 680f6dfbefdSBarry Smith Level: intermediate 6818f626970SPeter Brune 682f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 683f6dfbefdSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()` 6848f626970SPeter Brune @*/ 685d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESCompositeSetDamping(SNES snes, PetscInt n, PetscReal dmp) 686d71ae5a4SJacob Faibussowitsch { 6878f626970SPeter Brune PetscFunctionBegin; 6888f626970SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 689cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeSetDamping_C", (SNES, PetscInt, PetscReal), (snes, n, dmp)); 6903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6918f626970SPeter Brune } 6928f626970SPeter Brune 693d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_Composite(SNES snes) 694d71ae5a4SJacob Faibussowitsch { 6954ae48641SStefano Zampini Vec F, X, B, Y; 696eed5f15bSPeter Brune PetscInt i; 697b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 69872edecb9SPeter Brune SNESNormSchedule normtype; 699eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite *)snes->data; 700eed5f15bSPeter Brune 701eed5f15bSPeter Brune PetscFunctionBegin; 702eed5f15bSPeter Brune X = snes->vec_sol; 703eed5f15bSPeter Brune F = snes->vec_func; 704eed5f15bSPeter Brune B = snes->vec_rhs; 7054ae48641SStefano Zampini Y = snes->vec_sol_update; 706eed5f15bSPeter Brune 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 708eed5f15bSPeter Brune snes->iter = 0; 709eed5f15bSPeter Brune snes->norm = 0.; 7100b762f1fSPatrick Farrell comp->innerFailures = 0; 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 712eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 7139566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normtype)); 71472edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 715eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 7169566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 717eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 718eed5f15bSPeter Brune 719a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 7209566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 721a6da83ecSPatrick Farrell } else { 7229566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 723a6da83ecSPatrick Farrell } 724422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 7259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 726eed5f15bSPeter Brune snes->iter = 0; 727eed5f15bSPeter Brune snes->norm = fnorm; 7289566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7299566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 730eed5f15bSPeter Brune 731eed5f15bSPeter Brune /* test convergence */ 732d76a863bSStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 7332d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, snes->norm)); 7343ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 735eed5f15bSPeter Brune } else { 7369566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7379566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7389566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, snes->norm)); 739eed5f15bSPeter Brune } 740eed5f15bSPeter Brune 7414ae48641SStefano Zampini for (i = 0; i < snes->max_its; i++) { 742eed5f15bSPeter Brune /* Call general purpose update function */ 743dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 744eed5f15bSPeter Brune 745b22975d2SPatrick Farrell /* Copy the state before modification by application of the composite solver; 746b22975d2SPatrick Farrell we will subtract the new state after application */ 7479566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 748b22975d2SPatrick Farrell 749eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 7509566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_Additive(snes, X, B, F, &fnorm)); 751eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 7529566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_Multiplicative(snes, X, B, F, &fnorm)); 75390a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 7549566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_AdditiveOptimal(snes, X, B, F, &fnorm)); 7556c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported SNESComposite type"); 756d5a53f06SPeter Brune if (snes->reason < 0) break; 757d5a53f06SPeter Brune 758b22975d2SPatrick Farrell /* Compute the solution update for convergence testing */ 7599566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1.0, X)); 760b22975d2SPatrick Farrell 761eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 7629566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 763b22975d2SPatrick Farrell 764a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 7659566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 7669566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 7679566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 7689566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 7699566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 770a6da83ecSPatrick Farrell } else { 7719566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 7729566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 7739566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 774b22975d2SPatrick Farrell 7759566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 7769566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 7779566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 778a6da83ecSPatrick Farrell } 779422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 780b22975d2SPatrick Farrell } else if (normtype == SNES_NORM_ALWAYS) { 7819566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 7829566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 7839566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 7849566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 785eed5f15bSPeter Brune } 786eed5f15bSPeter Brune /* Monitor convergence */ 7879566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 788eed5f15bSPeter Brune snes->iter = i + 1; 789eed5f15bSPeter Brune snes->norm = fnorm; 790c1e67a49SFande Kong snes->xnorm = xnorm; 791c1e67a49SFande Kong snes->ynorm = snorm; 7929566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7939566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 794eed5f15bSPeter Brune /* Test for convergence */ 7952d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, snorm, fnorm)); 7962d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 797eed5f15bSPeter Brune if (snes->reason) break; 798eed5f15bSPeter Brune } 7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 800eed5f15bSPeter Brune } 801eed5f15bSPeter Brune 802eed5f15bSPeter Brune /*MC 803eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 804eed5f15bSPeter Brune 805eed5f15bSPeter Brune Options Database Keys: 806eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 8072fe279fdSBarry Smith - -snes_composite_sneses - <snes0,snes1,...> list of `SNES` to compose 808eed5f15bSPeter Brune 809eed5f15bSPeter Brune Level: intermediate 810eed5f15bSPeter Brune 8114f02bc6aSBarry Smith References: 812606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8134f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8144f02bc6aSBarry Smith 815f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 816f6dfbefdSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()`, 817f6dfbefdSBarry Smith `SNESCompositeSetDamping()` 818eed5f15bSPeter Brune M*/ 819eed5f15bSPeter Brune 820d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 821d71ae5a4SJacob Faibussowitsch { 822eed5f15bSPeter Brune SNES_Composite *jac; 823eed5f15bSPeter Brune 824eed5f15bSPeter Brune PetscFunctionBegin; 8254dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 826eed5f15bSPeter Brune 827eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 828eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 829eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 830eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 831eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 832eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 833eed5f15bSPeter Brune 834d8d34be6SBarry Smith snes->usesksp = PETSC_FALSE; 835d8d34be6SBarry Smith 8364fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8374fc747eaSLawrence Mitchell 838eed5f15bSPeter Brune snes->data = (void *)jac; 83990a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 84090a8ba9bSPeter Brune jac->Fes = NULL; 84190a8ba9bSPeter Brune jac->Xes = NULL; 8429c2f3473SPeter Brune jac->fnorms = NULL; 84390a8ba9bSPeter Brune jac->nsnes = 0; 8449e5d0892SLisandro Dalcin jac->head = NULL; 8455e806d2eSPeter Brune jac->stol = 0.1; 8465e806d2eSPeter Brune jac->rtol = 1.1; 847eed5f15bSPeter Brune 84890a8ba9bSPeter Brune jac->h = NULL; 84990a8ba9bSPeter Brune jac->s = NULL; 85090a8ba9bSPeter Brune jac->beta = NULL; 85190a8ba9bSPeter Brune jac->work = NULL; 85290a8ba9bSPeter Brune jac->rwork = NULL; 85390a8ba9bSPeter Brune 8549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", SNESCompositeSetType_Composite)); 8559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", SNESCompositeAddSNES_Composite)); 8569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", SNESCompositeGetSNES_Composite)); 8579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", SNESCompositeSetDamping_Composite)); 8583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 859eed5f15bSPeter Brune } 860