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 489371c9d4SSatish Balay static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) { 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; 63d5a53f06SPeter Brune PetscFunctionReturn(0); 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; 82d5a53f06SPeter Brune PetscFunctionReturn(0); 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 } 108eed5f15bSPeter Brune PetscFunctionReturn(0); 109eed5f15bSPeter Brune } 110eed5f15bSPeter Brune 1119371c9d4SSatish Balay static PetscErrorCode SNESCompositeApply_Additive(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) { 112eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 113eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 114eed5f15bSPeter Brune Vec Y, Xorig; 115d5a53f06SPeter Brune SNESConvergedReason reason; 116eed5f15bSPeter Brune 117eed5f15bSPeter Brune PetscFunctionBegin; 118eed5f15bSPeter Brune Y = snes->vec_sol_update; 1199566063dSJacob Faibussowitsch if (!jac->Xorig) PetscCall(VecDuplicate(X, &jac->Xorig)); 120eed5f15bSPeter Brune Xorig = jac->Xorig; 1219566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xorig)); 12228b400f6SJacob Faibussowitsch PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 12372edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 1249566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 125eed5f15bSPeter Brune while (next->next) { 126eed5f15bSPeter Brune next = next->next; 1279566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 128eed5f15bSPeter Brune } 129eed5f15bSPeter Brune } 130eed5f15bSPeter Brune next = jac->head; 1319566063dSJacob Faibussowitsch PetscCall(VecCopy(Xorig, Y)); 1329566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Y)); 1339566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 134d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1350b762f1fSPatrick Farrell jac->innerFailures++; 1360b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 137d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 138d5a53f06SPeter Brune PetscFunctionReturn(0); 139d5a53f06SPeter Brune } 1400b762f1fSPatrick Farrell } 1419566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, -1.0, Xorig)); 1429566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, next->dmp, Y)); 1438f626970SPeter Brune while (next->next) { 1448f626970SPeter Brune next = next->next; 1459566063dSJacob Faibussowitsch PetscCall(VecCopy(Xorig, Y)); 1469566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Y)); 1479566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 148d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1490b762f1fSPatrick Farrell jac->innerFailures++; 1500b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 151d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 152d5a53f06SPeter Brune PetscFunctionReturn(0); 153d5a53f06SPeter Brune } 1540b762f1fSPatrick Farrell } 1559566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, -1.0, Xorig)); 1569566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, next->dmp, Y)); 157eed5f15bSPeter Brune } 15872edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 1599566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 160a6da83ecSPatrick Farrell if (fnorm) { 161a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 1629566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 163a6da83ecSPatrick Farrell } else { 1649566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 165a6da83ecSPatrick Farrell } 166422a814eSBarry Smith SNESCheckFunctionNorm(snes, *fnorm); 167a6da83ecSPatrick Farrell } 168eed5f15bSPeter Brune } 169eed5f15bSPeter Brune PetscFunctionReturn(0); 170eed5f15bSPeter Brune } 17190a8ba9bSPeter Brune 17290a8ba9bSPeter Brune /* approximately solve the overdetermined system: 17390a8ba9bSPeter Brune 17490a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 17590a8ba9bSPeter Brune \alpha_i = 1 17690a8ba9bSPeter Brune 17790a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 17890a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 17990a8ba9bSPeter Brune 18090a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 18190a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 18290a8ba9bSPeter Brune */ 1839371c9d4SSatish Balay static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) { 18490a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 18590a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 18690a8ba9bSPeter Brune Vec *Xes = jac->Xes, *Fes = jac->Fes; 18790a8ba9bSPeter Brune PetscInt i, j; 1885e806d2eSPeter Brune PetscScalar tot, total, ftf; 1895e806d2eSPeter Brune PetscReal min_fnorm; 1905e806d2eSPeter Brune PetscInt min_i; 191d5a53f06SPeter Brune SNESConvergedReason reason; 19290a8ba9bSPeter Brune 19390a8ba9bSPeter Brune PetscFunctionBegin; 19428b400f6SJacob Faibussowitsch PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 19558bdfa14SPeter Brune 19658bdfa14SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 19758bdfa14SPeter Brune next = jac->head; 1989566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 19958bdfa14SPeter Brune while (next->next) { 20058bdfa14SPeter Brune next = next->next; 2019566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next->snes, F)); 20258bdfa14SPeter Brune } 20358bdfa14SPeter Brune } 20458bdfa14SPeter Brune 20590a8ba9bSPeter Brune next = jac->head; 20690a8ba9bSPeter Brune i = 0; 2079566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xes[i])); 2089566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Xes[i])); 2099566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 210d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2110b762f1fSPatrick Farrell jac->innerFailures++; 2120b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 213d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 214d5a53f06SPeter Brune PetscFunctionReturn(0); 215d5a53f06SPeter Brune } 2160b762f1fSPatrick Farrell } 21790a8ba9bSPeter Brune while (next->next) { 21890a8ba9bSPeter Brune i++; 21990a8ba9bSPeter Brune next = next->next; 2209566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xes[i])); 2219566063dSJacob Faibussowitsch PetscCall(SNESSolve(next->snes, B, Xes[i])); 2229566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next->snes, &reason)); 223d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2240b762f1fSPatrick Farrell jac->innerFailures++; 2250b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 226d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 227d5a53f06SPeter Brune PetscFunctionReturn(0); 228d5a53f06SPeter Brune } 22990a8ba9bSPeter Brune } 2300b762f1fSPatrick Farrell } 23190a8ba9bSPeter Brune 23290a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 23390a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 23448a46eb9SPierre Jolivet for (j = 0; j < i + 1; j++) PetscCall(VecDotBegin(Fes[i], Fes[j], &jac->h[i + j * jac->n])); 2359566063dSJacob Faibussowitsch PetscCall(VecDotBegin(Fes[i], F, &jac->g[i])); 23690a8ba9bSPeter Brune } 2375e806d2eSPeter Brune 23890a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 23990a8ba9bSPeter Brune for (j = 0; j < i + 1; j++) { 2409566063dSJacob Faibussowitsch PetscCall(VecDotEnd(Fes[i], Fes[j], &jac->h[i + j * jac->n])); 2419c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j * jac->n])); 24290a8ba9bSPeter Brune } 2439566063dSJacob Faibussowitsch PetscCall(VecDotEnd(Fes[i], F, &jac->g[i])); 24490a8ba9bSPeter Brune } 24590a8ba9bSPeter Brune 2469c2f3473SPeter Brune ftf = (*fnorm) * (*fnorm); 24790a8ba9bSPeter Brune 24890a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 249ad540459SPierre Jolivet for (j = i + 1; j < jac->n; j++) jac->h[i + j * jac->n] = jac->h[j + i * jac->n]; 25090a8ba9bSPeter Brune } 25190a8ba9bSPeter Brune 25290a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 253ad540459SPierre 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; 2549c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 2559c2f3473SPeter Brune } 25690a8ba9bSPeter Brune 25790a8ba9bSPeter Brune jac->info = 0; 25890a8ba9bSPeter Brune jac->rcond = -1.; 2599566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 26090a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 261792fecdfSBarry 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)); 26290a8ba9bSPeter Brune #else 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->info)); 26490a8ba9bSPeter Brune #endif 2659566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 26608401ef6SPierre Jolivet PetscCheck(jac->info >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "Bad argument to GELSS"); 26708401ef6SPierre Jolivet PetscCheck(jac->info <= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD failed to converge"); 2689c2f3473SPeter Brune tot = 0.; 2695e806d2eSPeter Brune total = 0.; 27090a8ba9bSPeter Brune for (i = 0; i < jac->n; i++) { 27108401ef6SPierre Jolivet PetscCheck(!snes->errorifnotconverged || !PetscIsInfOrNanScalar(jac->beta[i]), PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD generated inconsistent output"); 27263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "%" PetscInt_FMT ": %g\n", i, (double)PetscRealPart(jac->beta[i]))); 2739c2f3473SPeter Brune tot += jac->beta[i]; 2745e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 27590a8ba9bSPeter Brune } 2769566063dSJacob Faibussowitsch PetscCall(VecScale(X, (1. - tot))); 2779566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, jac->n, jac->beta, Xes)); 2789566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 279a6da83ecSPatrick Farrell 280a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 2819566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 282a6da83ecSPatrick Farrell } else { 2839566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 284a6da83ecSPatrick Farrell } 28590a8ba9bSPeter Brune 2865e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 2875e806d2eSPeter Brune min_fnorm = jac->fnorms[0]; 2885e806d2eSPeter Brune min_i = 0; 2899c2f3473SPeter Brune for (i = 0; i < jac->n; i++) { 2905e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) { 2915e806d2eSPeter Brune min_fnorm = jac->fnorms[i]; 2925e806d2eSPeter Brune min_i = i; 2939c2f3473SPeter Brune } 2949c2f3473SPeter Brune } 2955e806d2eSPeter Brune 2965e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */ 2975e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm * jac->rtol < *fnorm) { 2989566063dSJacob Faibussowitsch PetscCall(VecCopy(jac->Xes[min_i], X)); 2999566063dSJacob Faibussowitsch PetscCall(VecCopy(jac->Fes[min_i], F)); 3005e806d2eSPeter Brune *fnorm = min_fnorm; 3015e806d2eSPeter Brune } 30290a8ba9bSPeter Brune PetscFunctionReturn(0); 30390a8ba9bSPeter Brune } 30490a8ba9bSPeter Brune 3059371c9d4SSatish Balay static PetscErrorCode SNESSetUp_Composite(SNES snes) { 306dd515a93SPeter Brune DM dm; 307eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 308eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 30990a8ba9bSPeter Brune PetscInt n = 0, i; 31090a8ba9bSPeter Brune Vec F; 311eed5f15bSPeter Brune 312eed5f15bSPeter Brune PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 31463cdc2bdSPatrick Farrell 31563cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 31663cdc2bdSPatrick Farrell /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */ 3179566063dSJacob Faibussowitsch if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3189566063dSJacob Faibussowitsch if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 319dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu); 32063cdc2bdSPatrick Farrell } 32163cdc2bdSPatrick Farrell 322eed5f15bSPeter Brune while (next) { 32390a8ba9bSPeter Brune n++; 3249566063dSJacob Faibussowitsch PetscCall(SNESSetDM(next->snes, dm)); 3259566063dSJacob Faibussowitsch PetscCall(SNESSetApplicationContext(next->snes, snes->user)); 32663cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 32763cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 3289566063dSJacob Faibussowitsch PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds)); 32963cdc2bdSPatrick Farrell } else { 3309566063dSJacob Faibussowitsch PetscCall(SNESVISetVariableBounds(next->snes, snes->xl, snes->xu)); 33163cdc2bdSPatrick Farrell } 33263cdc2bdSPatrick Farrell } 33363cdc2bdSPatrick Farrell 334eed5f15bSPeter Brune next = next->next; 335eed5f15bSPeter Brune } 33690a8ba9bSPeter Brune jac->nsnes = n; 3379566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 33890a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 3399566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(F, jac->nsnes, &jac->Xes)); 3409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &jac->Fes)); 3419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &jac->fnorms)); 34290a8ba9bSPeter Brune next = jac->head; 34390a8ba9bSPeter Brune i = 0; 34490a8ba9bSPeter Brune while (next) { 3459566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(next->snes, &F, NULL, NULL)); 34690a8ba9bSPeter Brune jac->Fes[i] = F; 3479566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)F)); 34890a8ba9bSPeter Brune next = next->next; 34990a8ba9bSPeter Brune i++; 35090a8ba9bSPeter Brune } 35190a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 35290a8ba9bSPeter Brune jac->nrhs = 1; 3539c2f3473SPeter Brune jac->lda = jac->nsnes; 35490a8ba9bSPeter Brune jac->ldb = jac->nsnes; 35590a8ba9bSPeter Brune jac->n = jac->nsnes; 35690a8ba9bSPeter Brune 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n * jac->n, &jac->h)); 3589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n, &jac->beta)); 3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n, &jac->s)); 3609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n, &jac->g)); 3619c2f3473SPeter Brune jac->lwork = 12 * jac->n; 362dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 3639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->rwork)); 36490a8ba9bSPeter Brune #endif 3659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->lwork, &jac->work)); 36690a8ba9bSPeter Brune } 36790a8ba9bSPeter Brune 368eed5f15bSPeter Brune PetscFunctionReturn(0); 369eed5f15bSPeter Brune } 370eed5f15bSPeter Brune 3719371c9d4SSatish Balay static PetscErrorCode SNESReset_Composite(SNES snes) { 372eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 373eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 374eed5f15bSPeter Brune 375eed5f15bSPeter Brune PetscFunctionBegin; 376eed5f15bSPeter Brune while (next) { 3779566063dSJacob Faibussowitsch PetscCall(SNESReset(next->snes)); 378eed5f15bSPeter Brune next = next->next; 379eed5f15bSPeter Brune } 3809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->Xorig)); 3819566063dSJacob Faibussowitsch if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Xes)); 3829566063dSJacob Faibussowitsch if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Fes)); 3839566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->fnorms)); 3849566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->h)); 3859566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->s)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->g)); 3879566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->beta)); 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->work)); 3899566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->rwork)); 390eed5f15bSPeter Brune PetscFunctionReturn(0); 391eed5f15bSPeter Brune } 392eed5f15bSPeter Brune 3939371c9d4SSatish Balay static PetscErrorCode SNESDestroy_Composite(SNES snes) { 394eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 395eed5f15bSPeter Brune SNES_CompositeLink next = jac->head, next_tmp; 396eed5f15bSPeter Brune 397eed5f15bSPeter Brune PetscFunctionBegin; 3989566063dSJacob Faibussowitsch PetscCall(SNESReset_Composite(snes)); 399eed5f15bSPeter Brune while (next) { 4009566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&next->snes)); 401eed5f15bSPeter Brune next_tmp = next; 402eed5f15bSPeter Brune next = next->next; 4039566063dSJacob Faibussowitsch PetscCall(PetscFree(next_tmp)); 404eed5f15bSPeter Brune } 4052e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", NULL)); 4062e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", NULL)); 4072e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", NULL)); 4082e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", NULL)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 410eed5f15bSPeter Brune PetscFunctionReturn(0); 411eed5f15bSPeter Brune } 412eed5f15bSPeter Brune 4139371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_Composite(SNES snes, PetscOptionItems *PetscOptionsObject) { 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 } 445eed5f15bSPeter Brune PetscFunctionReturn(0); 446eed5f15bSPeter Brune } 447eed5f15bSPeter Brune 4489371c9d4SSatish Balay static PetscErrorCode SNESView_Composite(SNES snes, PetscViewer viewer) { 449eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 450eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 451eed5f15bSPeter Brune PetscBool iascii; 452eed5f15bSPeter Brune 453eed5f15bSPeter Brune PetscFunctionBegin; 4549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 455eed5f15bSPeter Brune if (iascii) { 4569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type - %s\n", SNESCompositeTypes[jac->type])); 4579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SNESes on composite preconditioner follow\n")); 4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 459eed5f15bSPeter Brune } 4601baa6e33SBarry Smith if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 461eed5f15bSPeter Brune while (next) { 4629566063dSJacob Faibussowitsch PetscCall(SNESView(next->snes, viewer)); 463eed5f15bSPeter Brune next = next->next; 464eed5f15bSPeter Brune } 465eed5f15bSPeter Brune if (iascii) { 4669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 4679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 468eed5f15bSPeter Brune } 469eed5f15bSPeter Brune PetscFunctionReturn(0); 470eed5f15bSPeter Brune } 471eed5f15bSPeter Brune 4729371c9d4SSatish Balay static PetscErrorCode SNESCompositeSetType_Composite(SNES snes, SNESCompositeType type) { 473eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite *)snes->data; 474eed5f15bSPeter Brune 475eed5f15bSPeter Brune PetscFunctionBegin; 476eed5f15bSPeter Brune jac->type = type; 477eed5f15bSPeter Brune PetscFunctionReturn(0); 478eed5f15bSPeter Brune } 479eed5f15bSPeter Brune 4809371c9d4SSatish Balay static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes, SNESType type) { 481eed5f15bSPeter Brune SNES_Composite *jac; 482eed5f15bSPeter Brune SNES_CompositeLink next, ilink; 483eed5f15bSPeter Brune PetscInt cnt = 0; 484eed5f15bSPeter Brune const char *prefix; 485d726e3a5SJed Brown char newprefix[20]; 486eed5f15bSPeter Brune DM dm; 487eed5f15bSPeter Brune 488eed5f15bSPeter Brune PetscFunctionBegin; 489*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ilink)); 4909e5d0892SLisandro Dalcin ilink->next = NULL; 4919566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &ilink->snes)); 4929566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 4939566063dSJacob Faibussowitsch PetscCall(SNESSetDM(ilink->snes, dm)); 4949566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(ilink->snes, snes->abstol, snes->rtol, snes->stol, 1, snes->max_funcs)); 4959566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)ilink->snes)); 496eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data; 497eed5f15bSPeter Brune next = jac->head; 498eed5f15bSPeter Brune if (!next) { 499eed5f15bSPeter Brune jac->head = ilink; 500eed5f15bSPeter Brune ilink->previous = NULL; 501eed5f15bSPeter Brune } else { 502eed5f15bSPeter Brune cnt++; 503eed5f15bSPeter Brune while (next->next) { 504eed5f15bSPeter Brune next = next->next; 505eed5f15bSPeter Brune cnt++; 506eed5f15bSPeter Brune } 507eed5f15bSPeter Brune next->next = ilink; 508eed5f15bSPeter Brune ilink->previous = next; 509eed5f15bSPeter Brune } 5109566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 5119566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(ilink->snes, prefix)); 5129566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(newprefix, sizeof(newprefix), "sub_%d_", (int)cnt)); 5139566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(ilink->snes, newprefix)); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes, (PetscObject)snes, 1)); 5159566063dSJacob Faibussowitsch PetscCall(SNESSetType(ilink->snes, type)); 5169566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY)); 51763cdc2bdSPatrick Farrell 5188f626970SPeter Brune ilink->dmp = 1.0; 51990a8ba9bSPeter Brune jac->nsnes++; 520eed5f15bSPeter Brune PetscFunctionReturn(0); 521eed5f15bSPeter Brune } 522eed5f15bSPeter Brune 5239371c9d4SSatish Balay static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes, PetscInt n, SNES *subsnes) { 524eed5f15bSPeter Brune SNES_Composite *jac; 525eed5f15bSPeter Brune SNES_CompositeLink next; 526eed5f15bSPeter Brune PetscInt i; 527eed5f15bSPeter Brune 528eed5f15bSPeter Brune PetscFunctionBegin; 529eed5f15bSPeter Brune jac = (SNES_Composite *)snes->data; 530eed5f15bSPeter Brune next = jac->head; 531eed5f15bSPeter Brune for (i = 0; i < n; i++) { 53228b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 533eed5f15bSPeter Brune next = next->next; 534eed5f15bSPeter Brune } 535eed5f15bSPeter Brune *subsnes = next->snes; 536eed5f15bSPeter Brune PetscFunctionReturn(0); 537eed5f15bSPeter Brune } 538eed5f15bSPeter Brune 539e31ab4f9SPeter Brune /*@C 540eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 541eed5f15bSPeter Brune 542f6dfbefdSBarry Smith Logically Collective on snes 543eed5f15bSPeter Brune 544d8d19677SJose E. Roman Input Parameters: 545eed5f15bSPeter Brune + snes - the preconditioner context 546f6dfbefdSBarry Smith - type - `SNES_COMPOSITE_ADDITIVE` (default), `SNES_COMPOSITE_MULTIPLICATIVE` 547eed5f15bSPeter Brune 548eed5f15bSPeter Brune Options Database Key: 549eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 550eed5f15bSPeter Brune 551eed5f15bSPeter Brune Level: Developer 552eed5f15bSPeter Brune 553f6dfbefdSBarry Smith .seealso: `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCOMPOSITE` 554eed5f15bSPeter Brune @*/ 5559371c9d4SSatish Balay PetscErrorCode SNESCompositeSetType(SNES snes, SNESCompositeType type) { 556eed5f15bSPeter Brune PetscFunctionBegin; 557eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 558eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes, type, 2); 559cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeSetType_C", (SNES, SNESCompositeType), (snes, type)); 560eed5f15bSPeter Brune PetscFunctionReturn(0); 561eed5f15bSPeter Brune } 562eed5f15bSPeter Brune 563eed5f15bSPeter Brune /*@C 564f6dfbefdSBarry Smith SNESCompositeAddSNES - Adds another `SNES` to the `SNESCOMPOSITE` 565eed5f15bSPeter Brune 566f6dfbefdSBarry Smith Collective on snes 567eed5f15bSPeter Brune 568eed5f15bSPeter Brune Input Parameters: 569f6dfbefdSBarry Smith + snes - the snes context of type `SNESCOMPOSITE` 570f6dfbefdSBarry Smith - type - the type of the new solver 571eed5f15bSPeter Brune 572eed5f15bSPeter Brune Level: Developer 573eed5f15bSPeter Brune 574f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeGetSNES()` 575eed5f15bSPeter Brune @*/ 5769371c9d4SSatish Balay PetscErrorCode SNESCompositeAddSNES(SNES snes, SNESType type) { 577eed5f15bSPeter Brune PetscFunctionBegin; 578eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 579cac4c232SBarry Smith PetscTryMethod(snes, "SNESCompositeAddSNES_C", (SNES, SNESType), (snes, type)); 580eed5f15bSPeter Brune PetscFunctionReturn(0); 581eed5f15bSPeter Brune } 5824ae48641SStefano Zampini 583eed5f15bSPeter Brune /*@ 584f6dfbefdSBarry Smith SNESCompositeGetSNES - Gets one of the `SNES` objects in the composite `SNESCOMPOSITE` 585eed5f15bSPeter Brune 586eed5f15bSPeter Brune Not Collective 587eed5f15bSPeter Brune 588d8d19677SJose E. Roman Input Parameters: 589f6dfbefdSBarry Smith + snes - the snes context 590f6dfbefdSBarry Smith - n - the number of the composed snes requested 591eed5f15bSPeter Brune 592eed5f15bSPeter Brune Output Parameters: 593f6dfbefdSBarry Smith . subsnes - the `SNES` requested 594eed5f15bSPeter Brune 595eed5f15bSPeter Brune Level: Developer 596eed5f15bSPeter Brune 597f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()` 598eed5f15bSPeter Brune @*/ 5999371c9d4SSatish Balay PetscErrorCode SNESCompositeGetSNES(SNES snes, PetscInt n, SNES *subsnes) { 600eed5f15bSPeter Brune PetscFunctionBegin; 601eed5f15bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 602eed5f15bSPeter Brune PetscValidPointer(subsnes, 3); 603cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeGetSNES_C", (SNES, PetscInt, SNES *), (snes, n, subsnes)); 604eed5f15bSPeter Brune PetscFunctionReturn(0); 605eed5f15bSPeter Brune } 606eed5f15bSPeter Brune 6076b2b7f7bSPatrick Farrell /*@ 608f6dfbefdSBarry Smith SNESCompositeGetNumber - Get the number of subsolvers in the `SNESCOMPOSITE` 6096b2b7f7bSPatrick Farrell 610f6dfbefdSBarry Smith Logically Collective on snes 6116b2b7f7bSPatrick Farrell 6126b2b7f7bSPatrick Farrell Input Parameter: 613f6dfbefdSBarry Smith snes - the snes context 6146b2b7f7bSPatrick Farrell 6156b2b7f7bSPatrick Farrell Output Parameter: 6166b2b7f7bSPatrick Farrell n - the number of subsolvers 6176b2b7f7bSPatrick Farrell 6186b2b7f7bSPatrick Farrell Level: Developer 6196b2b7f7bSPatrick Farrell 620f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()` 6216b2b7f7bSPatrick Farrell @*/ 6229371c9d4SSatish Balay PetscErrorCode SNESCompositeGetNumber(SNES snes, PetscInt *n) { 6236b2b7f7bSPatrick Farrell SNES_Composite *jac; 6246b2b7f7bSPatrick Farrell SNES_CompositeLink next; 6256b2b7f7bSPatrick Farrell 6266b2b7f7bSPatrick Farrell PetscFunctionBegin; 6276b2b7f7bSPatrick Farrell jac = (SNES_Composite *)snes->data; 6286b2b7f7bSPatrick Farrell next = jac->head; 6296b2b7f7bSPatrick Farrell 6306b2b7f7bSPatrick Farrell *n = 0; 6316b2b7f7bSPatrick Farrell while (next) { 6326b2b7f7bSPatrick Farrell *n = *n + 1; 6336b2b7f7bSPatrick Farrell next = next->next; 6346b2b7f7bSPatrick Farrell } 6356b2b7f7bSPatrick Farrell PetscFunctionReturn(0); 6366b2b7f7bSPatrick Farrell } 6376b2b7f7bSPatrick Farrell 6389371c9d4SSatish Balay static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes, PetscInt n, PetscReal dmp) { 6398f626970SPeter Brune SNES_Composite *jac; 6408f626970SPeter Brune SNES_CompositeLink next; 6418f626970SPeter Brune PetscInt i; 6428f626970SPeter Brune 6438f626970SPeter Brune PetscFunctionBegin; 6448f626970SPeter Brune jac = (SNES_Composite *)snes->data; 6458f626970SPeter Brune next = jac->head; 6468f626970SPeter Brune for (i = 0; i < n; i++) { 64728b400f6SJacob Faibussowitsch PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 6488f626970SPeter Brune next = next->next; 6498f626970SPeter Brune } 6508f626970SPeter Brune next->dmp = dmp; 6518f626970SPeter Brune PetscFunctionReturn(0); 6528f626970SPeter Brune } 6538f626970SPeter Brune 6548f626970SPeter Brune /*@ 655f6dfbefdSBarry Smith SNESCompositeSetDamping - Sets the damping of a subsolver when using `SNES_COMPOSITE_ADDITIVE` `SNESCOMPOSITE` 6568f626970SPeter Brune 6578f626970SPeter Brune Not Collective 6588f626970SPeter Brune 659d8d19677SJose E. Roman Input Parameters: 660f6dfbefdSBarry Smith + snes - the snes context 6618f626970SPeter Brune . n - the number of the snes requested 6628f626970SPeter Brune - dmp - the damping 6638f626970SPeter Brune 664f6dfbefdSBarry Smith Level: intermediate 6658f626970SPeter Brune 666f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 667f6dfbefdSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()` 6688f626970SPeter Brune @*/ 6699371c9d4SSatish Balay PetscErrorCode SNESCompositeSetDamping(SNES snes, PetscInt n, PetscReal dmp) { 6708f626970SPeter Brune PetscFunctionBegin; 6718f626970SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 672cac4c232SBarry Smith PetscUseMethod(snes, "SNESCompositeSetDamping_C", (SNES, PetscInt, PetscReal), (snes, n, dmp)); 6738f626970SPeter Brune PetscFunctionReturn(0); 6748f626970SPeter Brune } 6758f626970SPeter Brune 6769371c9d4SSatish Balay static PetscErrorCode SNESSolve_Composite(SNES snes) { 6774ae48641SStefano Zampini Vec F, X, B, Y; 678eed5f15bSPeter Brune PetscInt i; 679b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 68072edecb9SPeter Brune SNESNormSchedule normtype; 681eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite *)snes->data; 682eed5f15bSPeter Brune 683eed5f15bSPeter Brune PetscFunctionBegin; 684eed5f15bSPeter Brune X = snes->vec_sol; 685eed5f15bSPeter Brune F = snes->vec_func; 686eed5f15bSPeter Brune B = snes->vec_rhs; 6874ae48641SStefano Zampini Y = snes->vec_sol_update; 688eed5f15bSPeter Brune 6899566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 690eed5f15bSPeter Brune snes->iter = 0; 691eed5f15bSPeter Brune snes->norm = 0.; 6920b762f1fSPatrick Farrell comp->innerFailures = 0; 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 694eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 6959566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normtype)); 69672edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 697eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 6989566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 699eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 700eed5f15bSPeter Brune 701a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 7029566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 703a6da83ecSPatrick Farrell } else { 7049566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 705a6da83ecSPatrick Farrell } 706422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 708eed5f15bSPeter Brune snes->iter = 0; 709eed5f15bSPeter Brune snes->norm = fnorm; 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7119566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7129566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, snes->norm)); 713eed5f15bSPeter Brune 714eed5f15bSPeter Brune /* test convergence */ 715dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 716eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 717eed5f15bSPeter Brune } else { 7189566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7199566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7209566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, snes->norm)); 721eed5f15bSPeter Brune } 722eed5f15bSPeter Brune 7234ae48641SStefano Zampini for (i = 0; i < snes->max_its; i++) { 724eed5f15bSPeter Brune /* Call general purpose update function */ 725dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 726eed5f15bSPeter Brune 727b22975d2SPatrick Farrell /* Copy the state before modification by application of the composite solver; 728b22975d2SPatrick Farrell we will subtract the new state after application */ 7299566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 730b22975d2SPatrick Farrell 731eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 7329566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_Additive(snes, X, B, F, &fnorm)); 733eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 7349566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_Multiplicative(snes, X, B, F, &fnorm)); 73590a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 7369566063dSJacob Faibussowitsch PetscCall(SNESCompositeApply_AdditiveOptimal(snes, X, B, F, &fnorm)); 7376c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported SNESComposite type"); 738d5a53f06SPeter Brune if (snes->reason < 0) break; 739d5a53f06SPeter Brune 740b22975d2SPatrick Farrell /* Compute the solution update for convergence testing */ 7419566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1.0, X)); 742b22975d2SPatrick Farrell 743eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 7449566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 745b22975d2SPatrick Farrell 746a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 7479566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 7489566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 7499566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 7509566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 7519566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 752a6da83ecSPatrick Farrell } else { 7539566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 7549566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 7559566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 756b22975d2SPatrick Farrell 7579566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 7589566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 7599566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 760a6da83ecSPatrick Farrell } 761422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 762b22975d2SPatrick Farrell } else if (normtype == SNES_NORM_ALWAYS) { 7639566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 7649566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 7659566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 7669566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 767eed5f15bSPeter Brune } 768eed5f15bSPeter Brune /* Monitor convergence */ 7699566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 770eed5f15bSPeter Brune snes->iter = i + 1; 771eed5f15bSPeter Brune snes->norm = fnorm; 772c1e67a49SFande Kong snes->xnorm = xnorm; 773c1e67a49SFande Kong snes->ynorm = snorm; 7749566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7759566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7769566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 777eed5f15bSPeter Brune /* Test for convergence */ 778dbbe0bcdSBarry Smith if (normtype == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, snorm, fnorm, &snes->reason, snes->cnvP); 779eed5f15bSPeter Brune if (snes->reason) break; 780eed5f15bSPeter Brune } 78172edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS) { 782eed5f15bSPeter Brune if (i == snes->max_its) { 78363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its)); 784eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 785eed5f15bSPeter Brune } 786eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 787eed5f15bSPeter Brune PetscFunctionReturn(0); 788eed5f15bSPeter Brune } 789eed5f15bSPeter Brune 790eed5f15bSPeter Brune /*MC 791eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 792eed5f15bSPeter Brune 793eed5f15bSPeter Brune Options Database Keys: 794eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 795eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 796eed5f15bSPeter Brune 797eed5f15bSPeter Brune Level: intermediate 798eed5f15bSPeter Brune 7994f02bc6aSBarry Smith References: 800606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8014f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8024f02bc6aSBarry Smith 803f6dfbefdSBarry Smith .seealso: `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 804f6dfbefdSBarry Smith `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()`, 805f6dfbefdSBarry Smith `SNESCompositeSetDamping()` 806eed5f15bSPeter Brune M*/ 807eed5f15bSPeter Brune 8089371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) { 809eed5f15bSPeter Brune SNES_Composite *jac; 810eed5f15bSPeter Brune 811eed5f15bSPeter Brune PetscFunctionBegin; 812*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 813eed5f15bSPeter Brune 814eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 815eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 816eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 817eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 818eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 819eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 820eed5f15bSPeter Brune 821d8d34be6SBarry Smith snes->usesksp = PETSC_FALSE; 822d8d34be6SBarry Smith 8234fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8244fc747eaSLawrence Mitchell 825eed5f15bSPeter Brune snes->data = (void *)jac; 82690a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 82790a8ba9bSPeter Brune jac->Fes = NULL; 82890a8ba9bSPeter Brune jac->Xes = NULL; 8299c2f3473SPeter Brune jac->fnorms = NULL; 83090a8ba9bSPeter Brune jac->nsnes = 0; 8319e5d0892SLisandro Dalcin jac->head = NULL; 8325e806d2eSPeter Brune jac->stol = 0.1; 8335e806d2eSPeter Brune jac->rtol = 1.1; 834eed5f15bSPeter Brune 83590a8ba9bSPeter Brune jac->h = NULL; 83690a8ba9bSPeter Brune jac->s = NULL; 83790a8ba9bSPeter Brune jac->beta = NULL; 83890a8ba9bSPeter Brune jac->work = NULL; 83990a8ba9bSPeter Brune jac->rwork = NULL; 84090a8ba9bSPeter Brune 8419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", SNESCompositeSetType_Composite)); 8429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", SNESCompositeAddSNES_Composite)); 8439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", SNESCompositeGetSNES_Composite)); 8449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", SNESCompositeSetDamping_Composite)); 845eed5f15bSPeter Brune PetscFunctionReturn(0); 846eed5f15bSPeter Brune } 847