1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdm.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "MatMultASPIN" 6 PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y) 7 { 8 PetscErrorCode ierr; 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; 21 22 PetscFunctionBegin; 23 ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr); 24 snes = (SNES)ctx; 25 ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr); 26 ierr = SNESGetFunction(npc,&W,NULL,NULL);CHKERRQ(ierr); 27 ierr = PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);CHKERRQ(ierr); 28 if (!match) { 29 ierr = PetscObjectGetComm((PetscObject)snes,&comm); 30 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz"); 31 } 32 ierr = SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);CHKERRQ(ierr); 33 ierr = SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);CHKERRQ(ierr); 34 35 ierr = VecSet(Y,0);CHKERRQ(ierr); 36 ierr = MatMult(npc->jacobian_pre,X,W);CHKERRQ(ierr); 37 38 for (i=0;i<n;i++) { 39 ierr = VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 40 } 41 for (i=0;i<n;i++) { 42 ierr = VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 43 ierr = VecSet(x[i],0.);CHKERRQ(ierr); 44 ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr); 45 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 46 ierr = PCApply(pc,b[i],x[i]);CHKERRQ(ierr); 47 ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 48 } 49 for (i=0;i<n;i++) { 50 ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 51 } 52 PetscFunctionReturn(0); 53 } 54 55 #undef __FUNCT__ 56 #define __FUNCT__ "SNESCreate_ASPIN" 57 /* -------------------------------------------------------------------------- */ 58 /*MC 59 SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton 60 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 66 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: 69 70 \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b)) 71 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}. 75 76 Level: intermediate 77 78 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetPC(), SNESGetPCSide() 79 80 M*/ 81 PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes) 82 { 83 PetscErrorCode ierr; 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; 92 93 PetscFunctionBegin; 94 /* set up the solver */ 95 ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr); 96 ierr = SNESSetPCSide(snes,PC_LEFT);CHKERRQ(ierr); 97 ierr = SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED);CHKERRQ(ierr); 98 ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr); 99 ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr); 100 ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr); 101 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 102 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 103 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 104 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 105 ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 106 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);CHKERRQ(ierr); 107 108 /* set up the shell matrix */ 109 ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 110 ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr); 111 ierr = MatCreateShell(comm,n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr); 112 ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr); 113 ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr); 114 115 ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr); 116 ierr = MatDestroy(&aspinmat);CHKERRQ(ierr); 117 118 PetscFunctionReturn(0); 119 } 120