xref: /petsc/src/snes/impls/nasm/aspin.c (revision 6a4a1270853b5b1995c94de98f44d28bec57470a)
1 #include <petsc/private/snesimpl.h> /*I   "petscsnes.h"   I*/
2 #include <petscdm.h>
3 
4 PetscErrorCode MatMultASPIN(Mat m, Vec X, Vec Y) {
5   void       *ctx;
6   SNES        snes;
7   PetscInt    n, i;
8   VecScatter *oscatter;
9   SNES       *subsnes;
10   PetscBool   match;
11   MPI_Comm    comm;
12   KSP         ksp;
13   Vec        *x, *b;
14   Vec         W;
15   SNES        npc;
16   Mat         subJ, subpJ;
17 
18   PetscFunctionBegin;
19   PetscCall(MatShellGetContext(m, &ctx));
20   snes = (SNES)ctx;
21   PetscCall(SNESGetNPC(snes, &npc));
22   PetscCall(SNESGetFunction(npc, &W, NULL, NULL));
23   PetscCall(PetscObjectTypeCompare((PetscObject)npc, SNESNASM, &match));
24   if (!match) {
25     PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
26     SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
27   }
28   PetscCall(SNESNASMGetSubdomains(npc, &n, &subsnes, NULL, &oscatter, NULL));
29   PetscCall(SNESNASMGetSubdomainVecs(npc, &n, &x, &b, NULL, NULL));
30 
31   PetscCall(VecSet(Y, 0));
32   PetscCall(MatMult(npc->jacobian_pre, X, W));
33 
34   for (i = 0; i < n; i++) PetscCall(VecScatterBegin(oscatter[i], W, b[i], INSERT_VALUES, SCATTER_FORWARD));
35   for (i = 0; i < n; i++) {
36     PetscCall(VecScatterEnd(oscatter[i], W, b[i], INSERT_VALUES, SCATTER_FORWARD));
37     PetscCall(VecSet(x[i], 0.));
38     PetscCall(SNESGetJacobian(subsnes[i], &subJ, &subpJ, NULL, NULL));
39     PetscCall(SNESGetKSP(subsnes[i], &ksp));
40     PetscCall(KSPSetOperators(ksp, subJ, subpJ));
41     PetscCall(KSPSolve(ksp, b[i], x[i]));
42     PetscCall(VecScatterBegin(oscatter[i], x[i], Y, ADD_VALUES, SCATTER_REVERSE));
43     PetscCall(VecScatterEnd(oscatter[i], x[i], Y, ADD_VALUES, SCATTER_REVERSE));
44   }
45   PetscFunctionReturn(0);
46 }
47 
48 static PetscErrorCode SNESDestroy_ASPIN(SNES snes) {
49   PetscFunctionBegin;
50   PetscCall(SNESDestroy(&snes->npc));
51   /* reset NEWTONLS and free the data */
52   PetscCall(SNESReset(snes));
53   PetscCall(PetscFree(snes->data));
54   PetscFunctionReturn(0);
55 }
56 
57 /*MC
58       SNESASPIN - Helper `SNES` type for Additive-Schwarz Preconditioned Inexact Newton
59 
60    Options Database Keys:
61 +  -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type `NASM`)
62 .  -npc_sub_snes_ - options prefix of the subdomain nonlinear solves
63 .  -npc_sub_ksp_ - options prefix of the subdomain Krylov solver
64 -  -npc_sub_pc_ - options prefix of the subdomain preconditioner
65 
66     Notes:
67     This solver transform the given nonlinear problem to a new form and then runs matrix-free Newton-Krylov with no
68     preconditioner on that transformed problem.
69 
70     This routine sets up an instance of `SNESNETWONLS` with nonlinear left preconditioning.  It differs from other
71     similar functionality in `SNES` as it creates a linear shell matrix that corresponds to the product
72 
73     \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
74 
75     which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
76     nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
77     factorizations are reused on each application of J_b^{-1}.
78 
79     The Krylov method used in this nonlinear solver is run with NO preconditioner, because the preconditioning is done
80     at the nonlinear level, but the Jacobian for the original function must be provided (or calculated via coloring and
81     finite differences automatically) in the Pmat location of `SNESSetJacobian()` because the action of the original Jacobian
82     is needed by the shell matrix used to apply the Jacobian of the nonlinear preconditioned problem (see above).
83     Note that since the Pmat is not used to construct a preconditioner it could be provided in a matrix-free form.
84     The code for this implementation is a bit confusing because the Amat of `SNESSetJacobian()` applies the Jacobian of the
85     nonlinearly preconditioned function Jacobian while the Pmat provides the Jacobian of the original user provided function.
86     Note that the original `SNES` and nonlinear preconditioner preconditioner (see `SNESGetNPC()`), in this case `SNESNASM`, share
87     the same Jacobian matrices. `SNESNASM` computes the needed Jacobian in `SNESNASMComputeFinalJacobian_Private()`.
88 
89    Level: intermediate
90 
91    References:
92 +  * - X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms",  SIAM J. Sci. Comput., 24, 2002.
93 -  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
94    SIAM Review, 57(4), 2015
95 
96 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNASM`, `SNESGetNPC()`, `SNESGetNPCSide()`
97 
98 M*/
99 PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes) {
100   SNES           npc;
101   KSP            ksp;
102   PC             pc;
103   Mat            aspinmat;
104   Vec            F;
105   PetscInt       n;
106   SNESLineSearch linesearch;
107 
108   PetscFunctionBegin;
109   /* set up the solver */
110   PetscCall(SNESSetType(snes, SNESNEWTONLS));
111   PetscCall(SNESSetNPCSide(snes, PC_LEFT));
112   PetscCall(SNESSetFunctionType(snes, SNES_FUNCTION_PRECONDITIONED));
113   PetscCall(SNESGetNPC(snes, &npc));
114   PetscCall(SNESSetType(npc, SNESNASM));
115   PetscCall(SNESNASMSetType(npc, PC_ASM_BASIC));
116   PetscCall(SNESNASMSetComputeFinalJacobian(npc, PETSC_TRUE));
117   PetscCall(SNESGetKSP(snes, &ksp));
118   PetscCall(KSPGetPC(ksp, &pc));
119   PetscCall(PCSetType(pc, PCNONE));
120   PetscCall(SNESGetLineSearch(snes, &linesearch));
121   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
122 
123   /* set up the shell matrix */
124   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
125   PetscCall(VecGetLocalSize(F, &n));
126   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)snes), n, n, PETSC_DECIDE, PETSC_DECIDE, snes, &aspinmat));
127   PetscCall(MatSetType(aspinmat, MATSHELL));
128   PetscCall(MatShellSetOperation(aspinmat, MATOP_MULT, (void (*)(void))MatMultASPIN));
129   PetscCall(SNESSetJacobian(snes, aspinmat, NULL, NULL, NULL));
130   PetscCall(MatDestroy(&aspinmat));
131 
132   snes->ops->destroy = SNESDestroy_ASPIN;
133 
134   PetscFunctionReturn(0);
135 }
136