1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2d728fb7dSPeter Brune #include <petscdm.h> 3d728fb7dSPeter Brune 49371c9d4SSatish Balay PetscErrorCode MatMultASPIN(Mat m, Vec X, Vec Y) { 5d728fb7dSPeter Brune void *ctx; 6d728fb7dSPeter Brune SNES snes; 7d728fb7dSPeter Brune PetscInt n, i; 8d728fb7dSPeter Brune VecScatter *oscatter; 9d728fb7dSPeter Brune SNES *subsnes; 10d728fb7dSPeter Brune PetscBool match; 11d728fb7dSPeter Brune MPI_Comm comm; 12d728fb7dSPeter Brune KSP ksp; 13d728fb7dSPeter Brune Vec *x, *b; 14ed3c11a9SPeter Brune Vec W; 15d728fb7dSPeter Brune SNES npc; 16724f369dSPeter Brune Mat subJ, subpJ; 17d728fb7dSPeter Brune 18d728fb7dSPeter Brune PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(m, &ctx)); 20d728fb7dSPeter Brune snes = (SNES)ctx; 219566063dSJacob Faibussowitsch PetscCall(SNESGetNPC(snes, &npc)); 229566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(npc, &W, NULL, NULL)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)npc, SNESNASM, &match)); 24d728fb7dSPeter Brune if (!match) { 259566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 26d728fb7dSPeter Brune SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz"); 27d728fb7dSPeter Brune } 289566063dSJacob Faibussowitsch PetscCall(SNESNASMGetSubdomains(npc, &n, &subsnes, NULL, &oscatter, NULL)); 299566063dSJacob Faibussowitsch PetscCall(SNESNASMGetSubdomainVecs(npc, &n, &x, &b, NULL, NULL)); 30ed3c11a9SPeter Brune 319566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0)); 329566063dSJacob Faibussowitsch PetscCall(MatMult(npc->jacobian_pre, X, W)); 33ed3c11a9SPeter Brune 34*48a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(VecScatterBegin(oscatter[i], W, b[i], INSERT_VALUES, SCATTER_FORWARD)); 35d728fb7dSPeter Brune for (i = 0; i < n; i++) { 369566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscatter[i], W, b[i], INSERT_VALUES, SCATTER_FORWARD)); 379566063dSJacob Faibussowitsch PetscCall(VecSet(x[i], 0.)); 389566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(subsnes[i], &subJ, &subpJ, NULL, NULL)); 399566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(subsnes[i], &ksp)); 409566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, subJ, subpJ)); 419566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, b[i], x[i])); 429566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscatter[i], x[i], Y, ADD_VALUES, SCATTER_REVERSE)); 439566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscatter[i], x[i], Y, ADD_VALUES, SCATTER_REVERSE)); 44d728fb7dSPeter Brune } 45d728fb7dSPeter Brune PetscFunctionReturn(0); 46d728fb7dSPeter Brune } 47d728fb7dSPeter Brune 489371c9d4SSatish Balay static PetscErrorCode SNESDestroy_ASPIN(SNES snes) { 496cf69a2aSStefano Zampini PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes->npc)); 51263b3d60SStefano Zampini /* reset NEWTONLS and free the data */ 529566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 539566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 546cf69a2aSStefano Zampini PetscFunctionReturn(0); 556cf69a2aSStefano Zampini } 566cf69a2aSStefano Zampini 572e483eabSPeter Brune /* -------------------------------------------------------------------------- */ 582e483eabSPeter Brune /*MC 592e483eabSPeter Brune SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton 602e483eabSPeter Brune 612e483eabSPeter Brune Options Database: 622e483eabSPeter Brune + -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM) 632e483eabSPeter Brune . -npc_sub_snes_ - options prefix of the subdomain nonlinear solves 642e483eabSPeter Brune . -npc_sub_ksp_ - options prefix of the subdomain Krylov solver 652e483eabSPeter Brune - -npc_sub_pc_ - options prefix of the subdomain preconditioner 662e483eabSPeter Brune 67d14a3f37SRichard Tran Mills Notes: 68d14a3f37SRichard Tran Mills This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other 69d14a3f37SRichard Tran Mills similar functionality in SNES as it creates a linear shell matrix that corresponds to the product 702e483eabSPeter Brune 712e483eabSPeter Brune \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b)) 722e483eabSPeter Brune 732e483eabSPeter Brune which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of 742e483eabSPeter Brune nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner 752e483eabSPeter Brune factorizations are reused on each application of J_b^{-1}. 762e483eabSPeter Brune 77951fe5abSBarry Smith The Krylov method used in this nonlinear solver is run with NO preconditioner, because the preconditioning is done 78951fe5abSBarry Smith at the nonlinear level, but the Jacobian for the original function must be provided (or calculated via coloring and 79951fe5abSBarry Smith finite differences automatically) in the Pmat location of SNESSetJacobian() because the action of the original Jacobian 80951fe5abSBarry Smith is needed by the shell matrix used to apply the Jacobian of the nonlinear preconditioned problem (see above). 81951fe5abSBarry Smith Note that since the Pmat is not used to construct a preconditioner it could be provided in a matrix-free form. 82951fe5abSBarry Smith The code for this implementation is a bit confusing because the Amat of SNESSetJacobian() applies the Jacobian of the 83951fe5abSBarry Smith nonlinearly preconditioned function Jacobian while the Pmat provides the Jacobian of the original user provided function. 84951fe5abSBarry Smith Note that the original SNES and nonlinear preconditioner preconditioner (see SNESGetNPC()), in this case NASM, share 85e0c02066SRichard Tran Mills the same Jacobian matrices. SNESNASM computes the needed Jacobian in SNESNASMComputeFinalJacobian_Private(). 86951fe5abSBarry Smith 8706dd6b0eSSatish Balay Level: intermediate 882e483eabSPeter Brune 894f02bc6aSBarry Smith References: 90606c0280SSatish Balay + * - X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms", SIAM J. Sci. Comput., 24, 2002. 91606c0280SSatish Balay - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 924f02bc6aSBarry Smith SIAM Review, 57(4), 2015 934f02bc6aSBarry Smith 94db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNASM`, `SNESGetNPC()`, `SNESGetNPCSide()` 952e483eabSPeter Brune 962e483eabSPeter Brune M*/ 979371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes) { 98d728fb7dSPeter Brune SNES npc; 99d728fb7dSPeter Brune KSP ksp; 100d728fb7dSPeter Brune PC pc; 101d728fb7dSPeter Brune Mat aspinmat; 102d728fb7dSPeter Brune Vec F; 103d728fb7dSPeter Brune PetscInt n; 104ed3c11a9SPeter Brune SNESLineSearch linesearch; 105d728fb7dSPeter Brune 106d728fb7dSPeter Brune PetscFunctionBegin; 107d728fb7dSPeter Brune /* set up the solver */ 1089566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes, SNESNEWTONLS)); 1099566063dSJacob Faibussowitsch PetscCall(SNESSetNPCSide(snes, PC_LEFT)); 1109566063dSJacob Faibussowitsch PetscCall(SNESSetFunctionType(snes, SNES_FUNCTION_PRECONDITIONED)); 1119566063dSJacob Faibussowitsch PetscCall(SNESGetNPC(snes, &npc)); 1129566063dSJacob Faibussowitsch PetscCall(SNESSetType(npc, SNESNASM)); 1139566063dSJacob Faibussowitsch PetscCall(SNESNASMSetType(npc, PC_ASM_BASIC)); 1149566063dSJacob Faibussowitsch PetscCall(SNESNASMSetComputeFinalJacobian(npc, PETSC_TRUE)); 1159566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 1169566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1179566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 1189566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 119*48a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 120d728fb7dSPeter Brune 121d728fb7dSPeter Brune /* set up the shell matrix */ 1229566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 1239566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(F, &n)); 1249566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)snes), n, n, PETSC_DECIDE, PETSC_DECIDE, snes, &aspinmat)); 1259566063dSJacob Faibussowitsch PetscCall(MatSetType(aspinmat, MATSHELL)); 1269566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(aspinmat, MATOP_MULT, (void (*)(void))MatMultASPIN)); 1279566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, aspinmat, NULL, NULL, NULL)); 1289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aspinmat)); 129d728fb7dSPeter Brune 1306cf69a2aSStefano Zampini snes->ops->destroy = SNESDestroy_ASPIN; 1316cf69a2aSStefano Zampini 132d728fb7dSPeter Brune PetscFunctionReturn(0); 133d728fb7dSPeter Brune } 134