xref: /petsc/src/snes/impls/nasm/aspin.c (revision 6cf69a2a66989311d3c9fa42cfb2f1881ecfd724)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2d728fb7dSPeter Brune #include <petscdm.h>
3d728fb7dSPeter Brune 
4d728fb7dSPeter Brune #undef __FUNCT__
5d728fb7dSPeter Brune #define __FUNCT__ "MatMultASPIN"
6d728fb7dSPeter Brune PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y)
7d728fb7dSPeter Brune {
8d728fb7dSPeter Brune   PetscErrorCode ierr;
9d728fb7dSPeter Brune   void           *ctx;
10d728fb7dSPeter Brune   SNES           snes;
11d728fb7dSPeter Brune   PetscInt       n,i;
12d728fb7dSPeter Brune   VecScatter     *oscatter;
13d728fb7dSPeter Brune   SNES           *subsnes;
14d728fb7dSPeter Brune   PetscBool      match;
15d728fb7dSPeter Brune   MPI_Comm       comm;
16d728fb7dSPeter Brune   KSP            ksp;
17d728fb7dSPeter Brune   Vec            *x,*b;
18ed3c11a9SPeter Brune   Vec            W;
19d728fb7dSPeter Brune   SNES           npc;
20724f369dSPeter Brune   Mat            subJ,subpJ;
21d728fb7dSPeter Brune 
22d728fb7dSPeter Brune   PetscFunctionBegin;
23d728fb7dSPeter Brune   ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr);
24d728fb7dSPeter Brune   snes = (SNES)ctx;
25be95d8f1SBarry Smith   ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr);
26ed3c11a9SPeter Brune   ierr = SNESGetFunction(npc,&W,NULL,NULL);CHKERRQ(ierr);
27d728fb7dSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);CHKERRQ(ierr);
28d728fb7dSPeter Brune   if (!match) {
29302440fdSBarry Smith     ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
30d728fb7dSPeter Brune     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
31d728fb7dSPeter Brune   }
32d728fb7dSPeter Brune   ierr = SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);CHKERRQ(ierr);
33d728fb7dSPeter Brune   ierr = SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);CHKERRQ(ierr);
34ed3c11a9SPeter Brune 
35ed3c11a9SPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
36ed3c11a9SPeter Brune   ierr = MatMult(npc->jacobian_pre,X,W);CHKERRQ(ierr);
37ed3c11a9SPeter Brune 
38d728fb7dSPeter Brune   for (i=0;i<n;i++) {
39ed3c11a9SPeter Brune     ierr = VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
40d728fb7dSPeter Brune   }
41d728fb7dSPeter Brune   for (i=0;i<n;i++) {
42ed3c11a9SPeter Brune     ierr = VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
43ed3c11a9SPeter Brune     ierr = VecSet(x[i],0.);CHKERRQ(ierr);
44724f369dSPeter Brune     ierr = SNESGetJacobian(subsnes[i],&subJ,&subpJ,NULL,NULL);CHKERRQ(ierr);
45d728fb7dSPeter Brune     ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr);
46724f369dSPeter Brune     ierr = KSPSetOperators(ksp,subJ,subpJ);CHKERRQ(ierr);
47724f369dSPeter Brune     ierr = KSPSolve(ksp,b[i],x[i]);CHKERRQ(ierr);
48d728fb7dSPeter Brune     ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
49ed3c11a9SPeter Brune   }
50ed3c11a9SPeter Brune   for (i=0;i<n;i++) {
51d728fb7dSPeter Brune     ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
52d728fb7dSPeter Brune   }
53d728fb7dSPeter Brune   PetscFunctionReturn(0);
54d728fb7dSPeter Brune }
55d728fb7dSPeter Brune 
56d728fb7dSPeter Brune #undef __FUNCT__
57*6cf69a2aSStefano Zampini #define __FUNCT__ "SNESDestroy_ASPIN"
58*6cf69a2aSStefano Zampini PetscErrorCode SNESDestroy_ASPIN(SNES snes)
59*6cf69a2aSStefano Zampini {
60*6cf69a2aSStefano Zampini   PetscErrorCode ierr;
61*6cf69a2aSStefano Zampini 
62*6cf69a2aSStefano Zampini   PetscFunctionBegin;
63*6cf69a2aSStefano Zampini   ierr = SNESDestroy(&snes->pc);CHKERRQ(ierr);
64*6cf69a2aSStefano Zampini   PetscFunctionReturn(0);
65*6cf69a2aSStefano Zampini }
66*6cf69a2aSStefano Zampini 
67*6cf69a2aSStefano Zampini #undef __FUNCT__
68d728fb7dSPeter Brune #define __FUNCT__ "SNESCreate_ASPIN"
692e483eabSPeter Brune /* -------------------------------------------------------------------------- */
702e483eabSPeter Brune /*MC
712e483eabSPeter Brune       SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton
722e483eabSPeter Brune 
732e483eabSPeter Brune    Options Database:
742e483eabSPeter Brune +  -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM)
752e483eabSPeter Brune .  -npc_sub_snes_ - options prefix of the subdomain nonlinear solves
762e483eabSPeter Brune .  -npc_sub_ksp_ - options prefix of the subdomain Krylov solver
772e483eabSPeter Brune -  -npc_sub_pc_ - options prefix of the subdomain preconditioner
782e483eabSPeter Brune 
792e483eabSPeter Brune     Notes: This routine sets up an instance of NETWONLS with nonlinear left preconditioning.  It differs from other
802e483eabSPeter Brune     similar functionality in SNES as it creates a linear shell matrix that corresponds to the product:
812e483eabSPeter Brune 
822e483eabSPeter Brune     \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
832e483eabSPeter Brune 
842e483eabSPeter Brune     which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
852e483eabSPeter Brune     nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
862e483eabSPeter Brune     factorizations are reused on each application of J_b^{-1}.
872e483eabSPeter Brune 
882e483eabSPeter Brune    Level: intermediate
892e483eabSPeter Brune 
90be95d8f1SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetNPC(), SNESGetNPCSide()
912e483eabSPeter Brune 
922e483eabSPeter Brune M*/
938cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes)
94d728fb7dSPeter Brune {
95d728fb7dSPeter Brune   PetscErrorCode ierr;
96d728fb7dSPeter Brune   SNES           npc;
97d728fb7dSPeter Brune   KSP            ksp;
98d728fb7dSPeter Brune   PC             pc;
99d728fb7dSPeter Brune   Mat            aspinmat;
100d728fb7dSPeter Brune   Vec            F;
101d728fb7dSPeter Brune   PetscInt       n;
102ed3c11a9SPeter Brune   SNESLineSearch linesearch;
103d728fb7dSPeter Brune 
104d728fb7dSPeter Brune   PetscFunctionBegin;
105d728fb7dSPeter Brune   /* set up the solver */
106d728fb7dSPeter Brune   ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr);
107be95d8f1SBarry Smith   ierr = SNESSetNPCSide(snes,PC_LEFT);CHKERRQ(ierr);
108ce5a860aSPeter Brune   ierr = SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED);CHKERRQ(ierr);
109be95d8f1SBarry Smith   ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr);
110d728fb7dSPeter Brune   ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr);
111d78e9beeSPeter Brune   ierr = SNESNASMSetType(npc,PC_ASM_BASIC);CHKERRQ(ierr);
112d728fb7dSPeter Brune   ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr);
113d728fb7dSPeter Brune   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
114d728fb7dSPeter Brune   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
115d728fb7dSPeter Brune   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1167601faf0SJed Brown   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
117ed3c11a9SPeter Brune   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);CHKERRQ(ierr);
118d728fb7dSPeter Brune 
119d728fb7dSPeter Brune   /* set up the shell matrix */
120d728fb7dSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
121d728fb7dSPeter Brune   ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr);
122*6cf69a2aSStefano Zampini   ierr = MatCreateShell(PetscObjectComm((PetscObject)snes),n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr);
123d728fb7dSPeter Brune   ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr);
124d728fb7dSPeter Brune   ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr);
125d728fb7dSPeter Brune   ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr);
126d728fb7dSPeter Brune   ierr = MatDestroy(&aspinmat);CHKERRQ(ierr);
127d728fb7dSPeter Brune 
128*6cf69a2aSStefano Zampini   snes->ops->destroy = SNESDestroy_ASPIN;
129*6cf69a2aSStefano Zampini 
130d728fb7dSPeter Brune   PetscFunctionReturn(0);
131d728fb7dSPeter Brune }
132