Actual source code: aspin.c

petsc-3.4.5 2014-06-29
  1: #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
  2: #include <petscdm.h>

  6: PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y)
  7: {
  9:   void           *ctx;
 10:   SNES           snes;
 11:   PetscInt       n,i;
 12:   VecScatter     *oscatter;
 13:   SNES           *subsnes;
 14:   PetscBool      match;
 15:   MPI_Comm       comm;
 16:   KSP            ksp;
 17:   PC             pc;
 18:   Vec            *x,*b;
 19:   Vec            W;
 20:   SNES           npc;

 23:   MatShellGetContext(m,&ctx);
 24:   snes = (SNES)ctx;
 25:   SNESGetPC(snes,&npc);
 26:   SNESGetFunction(npc,&W,NULL,NULL);
 27:   PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);
 28:   if (!match) {
 29:     PetscObjectGetComm((PetscObject)snes,&comm);
 30:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
 31:   }
 32:   SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);
 33:   SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);

 35:   VecSet(Y,0);
 36:   MatMult(npc->jacobian_pre,X,W);

 38:   for (i=0;i<n;i++) {
 39:     VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);
 40:   }
 41:   for (i=0;i<n;i++) {
 42:     VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);
 43:     VecSet(x[i],0.);
 44:     SNESGetKSP(subsnes[i],&ksp);
 45:     KSPGetPC(ksp,&pc);
 46:     PCApply(pc,b[i],x[i]);
 47:     VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);
 48:   }
 49:   for (i=0;i<n;i++) {
 50:     VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);
 51:   }
 52:   return(0);
 53: }

 57: /* -------------------------------------------------------------------------- */
 58: /*MC
 59:       SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton

 61:    Options Database:
 62: +  -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM)
 63: .  -npc_sub_snes_ - options prefix of the subdomain nonlinear solves
 64: .  -npc_sub_ksp_ - options prefix of the subdomain Krylov solver
 65: -  -npc_sub_pc_ - options prefix of the subdomain preconditioner

 67:     Notes: This routine sets up an instance of NETWONLS with nonlinear left preconditioning.  It differs from other
 68:     similar functionality in SNES as it creates a linear shell matrix that corresponds to the product:

 70:     \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))

 72:     which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
 73:     nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
 74:     factorizations are reused on each application of J_b^{-1}.

 76:    Level: intermediate

 78: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetPC(), SNESGetPCSide()

 80: M*/
 81: PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes)
 82: {
 84:   SNES           npc;
 85:   KSP            ksp;
 86:   PC             pc;
 87:   Mat            aspinmat;
 88:   MPI_Comm       comm;
 89:   Vec            F;
 90:   PetscInt       n;
 91:   SNESLineSearch linesearch;

 94:   /* set up the solver */
 95:   SNESSetType(snes,SNESNEWTONLS);
 96:   SNESSetPCSide(snes,PC_LEFT);
 97:   SNESGetPC(snes,&npc);
 98:   SNESSetType(npc,SNESNASM);
 99:   SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);
100:   SNESGetKSP(snes,&ksp);
101:   KSPGetPC(ksp,&pc);
102:   PCSetType(pc,PCNONE);
103:   PetscObjectGetComm((PetscObject)snes,&comm);
104:   SNESGetLineSearch(snes,&linesearch);
105:   SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);

107:   /* set up the shell matrix */
108:   SNESGetFunction(snes,&F,NULL,NULL);
109:   VecGetLocalSize(F,&n);
110:   MatCreateShell(comm,n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);
111:   MatSetType(aspinmat,MATSHELL);
112:   MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);

114:   SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);
115:   MatDestroy(&aspinmat);

117:   return(0);
118: }