xref: /petsc/src/snes/impls/nasm/aspin.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2d728fb7dSPeter Brune #include <petscdm.h>
3d728fb7dSPeter Brune 
4d728fb7dSPeter Brune PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y)
5d728fb7dSPeter Brune {
6d728fb7dSPeter Brune   PetscErrorCode ierr;
7d728fb7dSPeter Brune   void           *ctx;
8d728fb7dSPeter Brune   SNES           snes;
9d728fb7dSPeter Brune   PetscInt       n,i;
10d728fb7dSPeter Brune   VecScatter     *oscatter;
11d728fb7dSPeter Brune   SNES           *subsnes;
12d728fb7dSPeter Brune   PetscBool      match;
13d728fb7dSPeter Brune   MPI_Comm       comm;
14d728fb7dSPeter Brune   KSP            ksp;
15d728fb7dSPeter Brune   Vec            *x,*b;
16ed3c11a9SPeter Brune   Vec            W;
17d728fb7dSPeter Brune   SNES           npc;
18724f369dSPeter Brune   Mat            subJ,subpJ;
19d728fb7dSPeter Brune 
20d728fb7dSPeter Brune   PetscFunctionBegin;
21d728fb7dSPeter Brune   ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr);
22d728fb7dSPeter Brune   snes = (SNES)ctx;
23be95d8f1SBarry Smith   ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr);
24ed3c11a9SPeter Brune   ierr = SNESGetFunction(npc,&W,NULL,NULL);CHKERRQ(ierr);
25d728fb7dSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);CHKERRQ(ierr);
26d728fb7dSPeter Brune   if (!match) {
27302440fdSBarry Smith     ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
28d728fb7dSPeter Brune     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
29d728fb7dSPeter Brune   }
30d728fb7dSPeter Brune   ierr = SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);CHKERRQ(ierr);
31d728fb7dSPeter Brune   ierr = SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);CHKERRQ(ierr);
32ed3c11a9SPeter Brune 
33ed3c11a9SPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
34ed3c11a9SPeter Brune   ierr = MatMult(npc->jacobian_pre,X,W);CHKERRQ(ierr);
35ed3c11a9SPeter Brune 
36d728fb7dSPeter Brune   for (i=0;i<n;i++) {
37ed3c11a9SPeter Brune     ierr = VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
38d728fb7dSPeter Brune   }
39d728fb7dSPeter Brune   for (i=0;i<n;i++) {
40ed3c11a9SPeter Brune     ierr = VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
41ed3c11a9SPeter Brune     ierr = VecSet(x[i],0.);CHKERRQ(ierr);
42724f369dSPeter Brune     ierr = SNESGetJacobian(subsnes[i],&subJ,&subpJ,NULL,NULL);CHKERRQ(ierr);
43d728fb7dSPeter Brune     ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr);
44724f369dSPeter Brune     ierr = KSPSetOperators(ksp,subJ,subpJ);CHKERRQ(ierr);
45724f369dSPeter Brune     ierr = KSPSolve(ksp,b[i],x[i]);CHKERRQ(ierr);
46d728fb7dSPeter Brune     ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
47ed3c11a9SPeter Brune   }
48ed3c11a9SPeter Brune   for (i=0;i<n;i++) {
49d728fb7dSPeter Brune     ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
50d728fb7dSPeter Brune   }
51d728fb7dSPeter Brune   PetscFunctionReturn(0);
52d728fb7dSPeter Brune }
53d728fb7dSPeter Brune 
54e0877f53SBarry Smith static PetscErrorCode SNESDestroy_ASPIN(SNES snes)
556cf69a2aSStefano Zampini {
566cf69a2aSStefano Zampini   PetscErrorCode ierr;
576cf69a2aSStefano Zampini 
586cf69a2aSStefano Zampini   PetscFunctionBegin;
59*efd4aadfSBarry Smith   ierr = SNESDestroy(&snes->npc);CHKERRQ(ierr);
60263b3d60SStefano Zampini   /* reset NEWTONLS and free the data */
61263b3d60SStefano Zampini   ierr = SNESReset(snes);CHKERRQ(ierr);
62263b3d60SStefano Zampini   ierr = PetscFree(snes->data);CHKERRQ(ierr);
636cf69a2aSStefano Zampini   PetscFunctionReturn(0);
646cf69a2aSStefano Zampini }
656cf69a2aSStefano Zampini 
662e483eabSPeter Brune /* -------------------------------------------------------------------------- */
672e483eabSPeter Brune /*MC
682e483eabSPeter Brune       SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton
692e483eabSPeter Brune 
702e483eabSPeter Brune    Options Database:
712e483eabSPeter Brune +  -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM)
722e483eabSPeter Brune .  -npc_sub_snes_ - options prefix of the subdomain nonlinear solves
732e483eabSPeter Brune .  -npc_sub_ksp_ - options prefix of the subdomain Krylov solver
742e483eabSPeter Brune -  -npc_sub_pc_ - options prefix of the subdomain preconditioner
752e483eabSPeter Brune 
762e483eabSPeter Brune     Notes: This routine sets up an instance of NETWONLS with nonlinear left preconditioning.  It differs from other
772e483eabSPeter Brune     similar functionality in SNES as it creates a linear shell matrix that corresponds to the product:
782e483eabSPeter Brune 
792e483eabSPeter Brune     \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
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
832e483eabSPeter Brune     factorizations are reused on each application of J_b^{-1}.
842e483eabSPeter Brune 
852e483eabSPeter Brune    Level: intermediate
862e483eabSPeter Brune 
874f02bc6aSBarry Smith    References:
8896a0c994SBarry Smith +  1. - X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms",  SIAM J. Sci. Comput., 24, 2002.
8996a0c994SBarry Smith -  2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
904f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
914f02bc6aSBarry Smith 
92be95d8f1SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetNPC(), SNESGetNPCSide()
932e483eabSPeter Brune 
942e483eabSPeter Brune M*/
958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes)
96d728fb7dSPeter Brune {
97d728fb7dSPeter Brune   PetscErrorCode ierr;
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 */
108d728fb7dSPeter Brune   ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr);
109be95d8f1SBarry Smith   ierr = SNESSetNPCSide(snes,PC_LEFT);CHKERRQ(ierr);
110ce5a860aSPeter Brune   ierr = SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED);CHKERRQ(ierr);
111be95d8f1SBarry Smith   ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr);
112d728fb7dSPeter Brune   ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr);
113d78e9beeSPeter Brune   ierr = SNESNASMSetType(npc,PC_ASM_BASIC);CHKERRQ(ierr);
114d728fb7dSPeter Brune   ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr);
115d728fb7dSPeter Brune   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
116d728fb7dSPeter Brune   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
117d728fb7dSPeter Brune   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1187601faf0SJed Brown   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
119ed3c11a9SPeter Brune   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);CHKERRQ(ierr);
120d728fb7dSPeter Brune 
121d728fb7dSPeter Brune   /* set up the shell matrix */
122d728fb7dSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
123d728fb7dSPeter Brune   ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr);
1246cf69a2aSStefano Zampini   ierr = MatCreateShell(PetscObjectComm((PetscObject)snes),n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr);
125d728fb7dSPeter Brune   ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr);
126d728fb7dSPeter Brune   ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr);
127d728fb7dSPeter Brune   ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr);
128d728fb7dSPeter Brune   ierr = MatDestroy(&aspinmat);CHKERRQ(ierr);
129d728fb7dSPeter Brune 
1306cf69a2aSStefano Zampini   snes->ops->destroy = SNESDestroy_ASPIN;
1316cf69a2aSStefano Zampini 
132d728fb7dSPeter Brune   PetscFunctionReturn(0);
133d728fb7dSPeter Brune }
134