xref: /petsc/src/snes/impls/nasm/aspin.c (revision 2286efddd54511ab18e8e2adb1e023c4bf8f0b92)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I   "petscsnes.h"   I*/
2d728fb7dSPeter Brune #include <petscdm.h>
3d728fb7dSPeter Brune 
MatMultASPIN(Mat m,Vec X,Vec Y)466976f2fSJacob Faibussowitsch static PetscErrorCode MatMultASPIN(Mat m, Vec X, Vec Y)
5d71ae5a4SJacob Faibussowitsch {
6d728fb7dSPeter Brune   void       *ctx;
7d728fb7dSPeter Brune   SNES        snes;
8d728fb7dSPeter Brune   PetscInt    n, i;
9d728fb7dSPeter Brune   VecScatter *oscatter;
10d728fb7dSPeter Brune   SNES       *subsnes;
11d728fb7dSPeter Brune   PetscBool   match;
12d728fb7dSPeter Brune   MPI_Comm    comm;
13d728fb7dSPeter Brune   KSP         ksp;
14d728fb7dSPeter Brune   Vec        *x, *b;
15ed3c11a9SPeter Brune   Vec         W;
16d728fb7dSPeter Brune   SNES        npc;
17724f369dSPeter Brune   Mat         subJ, subpJ;
18d728fb7dSPeter Brune 
19d728fb7dSPeter Brune   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(m, &ctx));
21d728fb7dSPeter Brune   snes = (SNES)ctx;
229566063dSJacob Faibussowitsch   PetscCall(SNESGetNPC(snes, &npc));
239566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(npc, &W, NULL, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)npc, SNESNASM, &match));
25d728fb7dSPeter Brune   if (!match) {
269566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
27d728fb7dSPeter Brune     SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
28d728fb7dSPeter Brune   }
299566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetSubdomains(npc, &n, &subsnes, NULL, &oscatter, NULL));
309566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetSubdomainVecs(npc, &n, &x, &b, NULL, NULL));
31ed3c11a9SPeter Brune 
329566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0));
339566063dSJacob Faibussowitsch   PetscCall(MatMult(npc->jacobian_pre, X, W));
34ed3c11a9SPeter Brune 
3548a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(VecScatterBegin(oscatter[i], W, b[i], INSERT_VALUES, SCATTER_FORWARD));
36d728fb7dSPeter Brune   for (i = 0; i < n; i++) {
379566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscatter[i], W, b[i], INSERT_VALUES, SCATTER_FORWARD));
389566063dSJacob Faibussowitsch     PetscCall(VecSet(x[i], 0.));
399566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(subsnes[i], &subJ, &subpJ, NULL, NULL));
409566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(subsnes[i], &ksp));
419566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, subJ, subpJ));
429566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ksp, b[i], x[i]));
439566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscatter[i], x[i], Y, ADD_VALUES, SCATTER_REVERSE));
449566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscatter[i], x[i], Y, ADD_VALUES, SCATTER_REVERSE));
45d728fb7dSPeter Brune   }
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47d728fb7dSPeter Brune }
48d728fb7dSPeter Brune 
SNESDestroy_ASPIN(SNES snes)49d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_ASPIN(SNES snes)
50d71ae5a4SJacob Faibussowitsch {
516cf69a2aSStefano Zampini   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes->npc));
53263b3d60SStefano Zampini   /* reset NEWTONLS and free the data */
549566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
576cf69a2aSStefano Zampini }
586cf69a2aSStefano Zampini 
592e483eabSPeter Brune /*MC
601d27aa22SBarry Smith    SNESASPIN - Helper `SNES` type for Additive-Schwarz Preconditioned Inexact Newton {cite}`ck02`, {cite}`bruneknepleysmithtu15`
612e483eabSPeter Brune 
62f6dfbefdSBarry Smith    Options Database Keys:
63f6dfbefdSBarry Smith +  -npc_snes_     - options prefix of the nonlinear subdomain solver (must be of type `NASM`)
642e483eabSPeter Brune .  -npc_sub_snes_ - options prefix of the subdomain nonlinear solves
652e483eabSPeter Brune .  -npc_sub_ksp_  - options prefix of the subdomain Krylov solver
662e483eabSPeter Brune -  -npc_sub_pc_   - options prefix of the subdomain preconditioner
672e483eabSPeter Brune 
68420bcc1bSBarry Smith      Level: intermediate
69420bcc1bSBarry Smith 
70d14a3f37SRichard Tran Mills     Notes:
71f6dfbefdSBarry Smith     This solver transform the given nonlinear problem to a new form and then runs matrix-free Newton-Krylov with no
72f6dfbefdSBarry Smith     preconditioner on that transformed problem.
73f6dfbefdSBarry Smith 
74f6dfbefdSBarry Smith     This routine sets up an instance of `SNESNETWONLS` with nonlinear left preconditioning.  It differs from other
75f6dfbefdSBarry Smith     similar functionality in `SNES` as it creates a linear shell matrix that corresponds to the product
762e483eabSPeter Brune 
771d27aa22SBarry Smith     $$
782e483eabSPeter Brune     \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
791d27aa22SBarry Smith     $$
802e483eabSPeter Brune 
812e483eabSPeter Brune     which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
822e483eabSPeter Brune     nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
831d27aa22SBarry Smith     factorizations are reused on each application of $J_b^{-1}$.
842e483eabSPeter Brune 
85951fe5abSBarry Smith     The Krylov method used in this nonlinear solver is run with NO preconditioner, because the preconditioning is done
86951fe5abSBarry Smith     at the nonlinear level, but the Jacobian for the original function must be provided (or calculated via coloring and
87f6dfbefdSBarry Smith     finite differences automatically) in the Pmat location of `SNESSetJacobian()` because the action of the original Jacobian
88951fe5abSBarry Smith     is needed by the shell matrix used to apply the Jacobian of the nonlinear preconditioned problem (see above).
89951fe5abSBarry Smith     Note that since the Pmat is not used to construct a preconditioner it could be provided in a matrix-free form.
90f6dfbefdSBarry Smith     The code for this implementation is a bit confusing because the Amat of `SNESSetJacobian()` applies the Jacobian of the
91951fe5abSBarry Smith     nonlinearly preconditioned function Jacobian while the Pmat provides the Jacobian of the original user provided function.
9215229ffcSPierre Jolivet     Note that the original `SNES` and nonlinear preconditioner (see `SNESGetNPC()`), in this case `SNESNASM`, share
93f6dfbefdSBarry Smith     the same Jacobian matrices. `SNESNASM` computes the needed Jacobian in `SNESNASMComputeFinalJacobian_Private()`.
94951fe5abSBarry Smith 
95420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNASM`, `SNESGetNPC()`, `SNESGetNPCSide()`
962e483eabSPeter Brune M*/
SNESCreate_ASPIN(SNES snes)97d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes)
98d71ae5a4SJacob Faibussowitsch {
99d728fb7dSPeter Brune   SNES           npc;
100d728fb7dSPeter Brune   KSP            ksp;
101d728fb7dSPeter Brune   PC             pc;
102d728fb7dSPeter Brune   Mat            aspinmat;
103d728fb7dSPeter Brune   Vec            F;
104d728fb7dSPeter Brune   PetscInt       n;
105ed3c11a9SPeter Brune   SNESLineSearch linesearch;
106d728fb7dSPeter Brune 
107d728fb7dSPeter Brune   PetscFunctionBegin;
108d728fb7dSPeter Brune   /* set up the solver */
1099566063dSJacob Faibussowitsch   PetscCall(SNESSetType(snes, SNESNEWTONLS));
1109566063dSJacob Faibussowitsch   PetscCall(SNESSetNPCSide(snes, PC_LEFT));
1119566063dSJacob Faibussowitsch   PetscCall(SNESSetFunctionType(snes, SNES_FUNCTION_PRECONDITIONED));
1129566063dSJacob Faibussowitsch   PetscCall(SNESGetNPC(snes, &npc));
1139566063dSJacob Faibussowitsch   PetscCall(SNESSetType(npc, SNESNASM));
1149566063dSJacob Faibussowitsch   PetscCall(SNESNASMSetType(npc, PC_ASM_BASIC));
1159566063dSJacob Faibussowitsch   PetscCall(SNESNASMSetComputeFinalJacobian(npc, PETSC_TRUE));
1169566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
1179566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
1189566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
1199566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
12048a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
121d728fb7dSPeter Brune 
122d728fb7dSPeter Brune   /* set up the shell matrix */
1239566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
1249566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(F, &n));
1259566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)snes), n, n, PETSC_DECIDE, PETSC_DECIDE, snes, &aspinmat));
1269566063dSJacob Faibussowitsch   PetscCall(MatSetType(aspinmat, MATSHELL));
127*57d50842SBarry Smith   PetscCall(MatShellSetOperation(aspinmat, MATOP_MULT, (PetscErrorCodeFn *)MatMultASPIN));
1289566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, aspinmat, NULL, NULL, NULL));
1299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aspinmat));
130d728fb7dSPeter Brune 
1316cf69a2aSStefano Zampini   snes->ops->destroy = SNESDestroy_ASPIN;
1322891b2c6SStefano Zampini   PetscCall(PetscObjectChangeTypeName((PetscObject)snes, SNESASPIN));
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134d728fb7dSPeter Brune }
135