xref: /petsc/src/snes/impls/nasm/aspin.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
1d728fb7dSPeter Brune #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;
17ca44f815SPeter Brune   PC             pc;
18d728fb7dSPeter Brune   Vec            *x,*b;
19ed3c11a9SPeter Brune   Vec            W;
20d728fb7dSPeter Brune   SNES           npc;
21d728fb7dSPeter Brune 
22d728fb7dSPeter Brune   PetscFunctionBegin;
23d728fb7dSPeter Brune   ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr);
24d728fb7dSPeter Brune   snes = (SNES)ctx;
25d728fb7dSPeter Brune   ierr = SNESGetPC(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) {
29d728fb7dSPeter Brune     ierr = PetscObjectGetComm((PetscObject)snes,&comm);
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);
44d728fb7dSPeter Brune     ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr);
45ca44f815SPeter Brune     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
46ca44f815SPeter Brune     ierr = PCApply(pc,b[i],x[i]);CHKERRQ(ierr);
47d728fb7dSPeter Brune     ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
48ed3c11a9SPeter Brune   }
49ed3c11a9SPeter Brune   for (i=0;i<n;i++) {
50d728fb7dSPeter Brune     ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
51d728fb7dSPeter Brune   }
52d728fb7dSPeter Brune   PetscFunctionReturn(0);
53d728fb7dSPeter Brune }
54d728fb7dSPeter Brune 
55d728fb7dSPeter Brune #undef __FUNCT__
56d728fb7dSPeter Brune #define __FUNCT__ "SNESCreate_ASPIN"
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 
672e483eabSPeter Brune     Notes: This routine sets up an instance of NETWONLS with nonlinear left preconditioning.  It differs from other
682e483eabSPeter Brune     similar functionality in SNES as it creates a linear shell matrix that corresponds to the product:
692e483eabSPeter Brune 
702e483eabSPeter Brune     \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
712e483eabSPeter Brune 
722e483eabSPeter Brune     which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
732e483eabSPeter Brune     nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
742e483eabSPeter Brune     factorizations are reused on each application of J_b^{-1}.
752e483eabSPeter Brune 
762e483eabSPeter Brune    Level: intermediate
772e483eabSPeter Brune 
782e483eabSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetPC(), SNESGetPCSide()
792e483eabSPeter Brune 
802e483eabSPeter Brune M*/
81*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes)
82d728fb7dSPeter Brune {
83d728fb7dSPeter Brune   PetscErrorCode ierr;
84d728fb7dSPeter Brune   SNES           npc;
85d728fb7dSPeter Brune   KSP            ksp;
86d728fb7dSPeter Brune   PC             pc;
87d728fb7dSPeter Brune   Mat            aspinmat;
88d728fb7dSPeter Brune   MPI_Comm       comm;
89d728fb7dSPeter Brune   Vec            F;
90d728fb7dSPeter Brune   PetscInt       n;
91ed3c11a9SPeter Brune   SNESLineSearch linesearch;
92d728fb7dSPeter Brune 
93d728fb7dSPeter Brune   PetscFunctionBegin;
94d728fb7dSPeter Brune   /* set up the solver */
95d728fb7dSPeter Brune   ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr);
96d728fb7dSPeter Brune   ierr = SNESSetPCSide(snes,PC_LEFT);CHKERRQ(ierr);
97d728fb7dSPeter Brune   ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr);
98d728fb7dSPeter Brune   ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr);
99d728fb7dSPeter Brune   ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr);
100d728fb7dSPeter Brune   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
101d728fb7dSPeter Brune   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
102d728fb7dSPeter Brune   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
103ed3c11a9SPeter Brune   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
104ed3c11a9SPeter Brune   ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
105ed3c11a9SPeter Brune   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);CHKERRQ(ierr);
106d728fb7dSPeter Brune 
107d728fb7dSPeter Brune   /* set up the shell matrix */
108d728fb7dSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
109d728fb7dSPeter Brune   ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr);
110d728fb7dSPeter Brune   ierr = MatCreateShell(comm,n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr);
111d728fb7dSPeter Brune   ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr);
112d728fb7dSPeter Brune   ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr);
113d728fb7dSPeter Brune 
114d728fb7dSPeter Brune   ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr);
115d728fb7dSPeter Brune   ierr = MatDestroy(&aspinmat);CHKERRQ(ierr);
116d728fb7dSPeter Brune 
117d728fb7dSPeter Brune   PetscFunctionReturn(0);
118d728fb7dSPeter Brune }
119