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; 17d728fb7dSPeter Brune Vec *x,*b; 18d728fb7dSPeter Brune SNES npc; 19d728fb7dSPeter Brune 20d728fb7dSPeter Brune PetscFunctionBegin; 21d728fb7dSPeter Brune ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr); 22d728fb7dSPeter Brune snes = (SNES)ctx; 23d728fb7dSPeter Brune ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr); 24d728fb7dSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);CHKERRQ(ierr); 25d728fb7dSPeter Brune if (!match) { 26d728fb7dSPeter Brune ierr = PetscObjectGetComm((PetscObject)snes,&comm); 27d728fb7dSPeter Brune SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz"); 28d728fb7dSPeter Brune } 29d728fb7dSPeter Brune ierr = SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);CHKERRQ(ierr); 30d728fb7dSPeter Brune ierr = SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);CHKERRQ(ierr); 31d728fb7dSPeter Brune ierr = MatMult(npc->jacobian_pre,X,Y);CHKERRQ(ierr); 32d728fb7dSPeter Brune for (i=0;i<n;i++) { 33d728fb7dSPeter Brune ierr = VecScatterBegin(oscatter[i],Y,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34d728fb7dSPeter Brune ierr = VecScatterEnd(oscatter[i],Y,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35d728fb7dSPeter Brune } 36d728fb7dSPeter Brune ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 37d728fb7dSPeter Brune for (i=0;i<n;i++) { 38d728fb7dSPeter Brune ierr = VecSet(x[i],0);CHKERRQ(ierr); 39d728fb7dSPeter Brune ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr); 40d728fb7dSPeter Brune ierr = KSPSolve(ksp,b[i],x[i]);CHKERRQ(ierr); 41d728fb7dSPeter Brune } 42d728fb7dSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 43d728fb7dSPeter Brune ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 44d728fb7dSPeter Brune for (i=0;i<n;i++) { 45d728fb7dSPeter Brune ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 46d728fb7dSPeter Brune ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 47d728fb7dSPeter Brune } 48d728fb7dSPeter Brune 49d728fb7dSPeter Brune ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 50d728fb7dSPeter Brune 51d728fb7dSPeter Brune PetscFunctionReturn(0); 52d728fb7dSPeter Brune } 53d728fb7dSPeter Brune 54d728fb7dSPeter Brune EXTERN_C_BEGIN 55d728fb7dSPeter Brune #undef __FUNCT__ 56d728fb7dSPeter Brune #define __FUNCT__ "SNESCreate_ASPIN" 57*2e483eabSPeter Brune /* -------------------------------------------------------------------------- */ 58*2e483eabSPeter Brune /*MC 59*2e483eabSPeter Brune SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton 60*2e483eabSPeter Brune 61*2e483eabSPeter Brune Options Database: 62*2e483eabSPeter Brune + -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM) 63*2e483eabSPeter Brune . -npc_sub_snes_ - options prefix of the subdomain nonlinear solves 64*2e483eabSPeter Brune . -npc_sub_ksp_ - options prefix of the subdomain Krylov solver 65*2e483eabSPeter Brune - -npc_sub_pc_ - options prefix of the subdomain preconditioner 66*2e483eabSPeter Brune 67*2e483eabSPeter Brune Notes: This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other 68*2e483eabSPeter Brune similar functionality in SNES as it creates a linear shell matrix that corresponds to the product: 69*2e483eabSPeter Brune 70*2e483eabSPeter Brune \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b)) 71*2e483eabSPeter Brune 72*2e483eabSPeter Brune which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of 73*2e483eabSPeter Brune nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner 74*2e483eabSPeter Brune factorizations are reused on each application of J_b^{-1}. 75*2e483eabSPeter Brune 76*2e483eabSPeter Brune Level: intermediate 77*2e483eabSPeter Brune 78*2e483eabSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetPC(), SNESGetPCSide() 79*2e483eabSPeter Brune 80*2e483eabSPeter Brune M*/ 81d728fb7dSPeter Brune 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; 91d728fb7dSPeter Brune 92d728fb7dSPeter Brune PetscFunctionBegin; 93d728fb7dSPeter Brune /* set up the solver */ 94d728fb7dSPeter Brune ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr); 95d728fb7dSPeter Brune ierr = SNESSetPCSide(snes,PC_LEFT);CHKERRQ(ierr); 96d728fb7dSPeter Brune ierr = SNESGetPC(snes,&npc);CHKERRQ(ierr); 97d728fb7dSPeter Brune ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr); 98d728fb7dSPeter Brune ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr); 99d728fb7dSPeter Brune ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 100d728fb7dSPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 101d728fb7dSPeter Brune ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 102d728fb7dSPeter Brune ierr = PetscObjectGetComm((PetscObject)snes,&comm); 103d728fb7dSPeter Brune 104d728fb7dSPeter Brune /* set up the shell matrix */ 105d728fb7dSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 106d728fb7dSPeter Brune ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr); 107d728fb7dSPeter Brune ierr = MatCreateShell(comm,n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr); 108d728fb7dSPeter Brune ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr); 109d728fb7dSPeter Brune ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr); 110d728fb7dSPeter Brune 111d728fb7dSPeter Brune ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr); 112d728fb7dSPeter Brune ierr = MatDestroy(&aspinmat);CHKERRQ(ierr); 113d728fb7dSPeter Brune 114d728fb7dSPeter Brune PetscFunctionReturn(0); 115d728fb7dSPeter Brune } 116d728fb7dSPeter Brune EXTERN_C_END 117