xref: /petsc/src/snes/impls/nasm/aspin.c (revision d728fb7dff1ad61d29df6dcab893bdd59c1da584)
1*d728fb7dSPeter Brune #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2*d728fb7dSPeter Brune #include <petscdm.h>
3*d728fb7dSPeter Brune 
4*d728fb7dSPeter Brune #undef __FUNCT__
5*d728fb7dSPeter Brune #define __FUNCT__ "MatMultASPIN"
6*d728fb7dSPeter Brune PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y)
7*d728fb7dSPeter Brune {
8*d728fb7dSPeter Brune   PetscErrorCode ierr;
9*d728fb7dSPeter Brune   void           *ctx;
10*d728fb7dSPeter Brune   SNES           snes;
11*d728fb7dSPeter Brune   PetscInt       n,i;
12*d728fb7dSPeter Brune   VecScatter     *oscatter;
13*d728fb7dSPeter Brune   SNES           *subsnes;
14*d728fb7dSPeter Brune   PetscBool      match;
15*d728fb7dSPeter Brune   MPI_Comm       comm;
16*d728fb7dSPeter Brune   KSP            ksp;
17*d728fb7dSPeter Brune   Vec            *x,*b;
18*d728fb7dSPeter Brune   SNES           npc;
19*d728fb7dSPeter Brune 
20*d728fb7dSPeter Brune   PetscFunctionBegin;
21*d728fb7dSPeter Brune   ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr);
22*d728fb7dSPeter Brune   snes = (SNES)ctx;
23*d728fb7dSPeter Brune   ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr);
24*d728fb7dSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);CHKERRQ(ierr);
25*d728fb7dSPeter Brune   if (!match) {
26*d728fb7dSPeter Brune     ierr = PetscObjectGetComm((PetscObject)snes,&comm);
27*d728fb7dSPeter Brune     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
28*d728fb7dSPeter Brune   }
29*d728fb7dSPeter Brune   ierr = SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);CHKERRQ(ierr);
30*d728fb7dSPeter Brune   ierr = SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);CHKERRQ(ierr);
31*d728fb7dSPeter Brune   ierr = MatMult(npc->jacobian_pre,X,Y);CHKERRQ(ierr);
32*d728fb7dSPeter Brune   for (i=0;i<n;i++) {
33*d728fb7dSPeter Brune     ierr = VecScatterBegin(oscatter[i],Y,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34*d728fb7dSPeter Brune     ierr = VecScatterEnd(oscatter[i],Y,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
35*d728fb7dSPeter Brune   }
36*d728fb7dSPeter Brune   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
37*d728fb7dSPeter Brune   for (i=0;i<n;i++) {
38*d728fb7dSPeter Brune     ierr = VecSet(x[i],0);CHKERRQ(ierr);
39*d728fb7dSPeter Brune     ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr);
40*d728fb7dSPeter Brune     ierr = KSPSolve(ksp,b[i],x[i]);CHKERRQ(ierr);
41*d728fb7dSPeter Brune   }
42*d728fb7dSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
43*d728fb7dSPeter Brune   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
44*d728fb7dSPeter Brune   for (i=0;i<n;i++) {
45*d728fb7dSPeter Brune     ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
46*d728fb7dSPeter Brune     ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
47*d728fb7dSPeter Brune   }
48*d728fb7dSPeter Brune 
49*d728fb7dSPeter Brune   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
50*d728fb7dSPeter Brune 
51*d728fb7dSPeter Brune   PetscFunctionReturn(0);
52*d728fb7dSPeter Brune }
53*d728fb7dSPeter Brune 
54*d728fb7dSPeter Brune EXTERN_C_BEGIN
55*d728fb7dSPeter Brune #undef __FUNCT__
56*d728fb7dSPeter Brune #define __FUNCT__ "SNESCreate_ASPIN"
57*d728fb7dSPeter Brune PetscErrorCode SNESCreate_ASPIN(SNES snes)
58*d728fb7dSPeter Brune {
59*d728fb7dSPeter Brune   PetscErrorCode ierr;
60*d728fb7dSPeter Brune   SNES           npc;
61*d728fb7dSPeter Brune   KSP            ksp;
62*d728fb7dSPeter Brune   PC             pc;
63*d728fb7dSPeter Brune   Mat            aspinmat;
64*d728fb7dSPeter Brune   MPI_Comm       comm;
65*d728fb7dSPeter Brune   Vec            F;
66*d728fb7dSPeter Brune   PetscInt       n;
67*d728fb7dSPeter Brune 
68*d728fb7dSPeter Brune   PetscFunctionBegin;
69*d728fb7dSPeter Brune   /* set up the solver */
70*d728fb7dSPeter Brune   ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr);
71*d728fb7dSPeter Brune   ierr = SNESSetPCSide(snes,PC_LEFT);CHKERRQ(ierr);
72*d728fb7dSPeter Brune   ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr);
73*d728fb7dSPeter Brune   ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr);
74*d728fb7dSPeter Brune   ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr);
75*d728fb7dSPeter Brune   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
76*d728fb7dSPeter Brune   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
77*d728fb7dSPeter Brune   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
78*d728fb7dSPeter Brune   ierr = PetscObjectGetComm((PetscObject)snes,&comm);
79*d728fb7dSPeter Brune 
80*d728fb7dSPeter Brune   /* set up the shell matrix */
81*d728fb7dSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
82*d728fb7dSPeter Brune   ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr);
83*d728fb7dSPeter Brune   ierr = MatCreateShell(comm,n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr);
84*d728fb7dSPeter Brune   ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr);
85*d728fb7dSPeter Brune   ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr);
86*d728fb7dSPeter Brune 
87*d728fb7dSPeter Brune   ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr);
88*d728fb7dSPeter Brune   ierr = MatDestroy(&aspinmat);CHKERRQ(ierr);
89*d728fb7dSPeter Brune 
90*d728fb7dSPeter Brune   PetscFunctionReturn(0);
91*d728fb7dSPeter Brune }
92*d728fb7dSPeter Brune EXTERN_C_END
93