1af0996ceSBarry Smith #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; 18ed3c11a9SPeter Brune Vec W; 19d728fb7dSPeter Brune SNES npc; 20724f369dSPeter Brune Mat subJ,subpJ; 21d728fb7dSPeter Brune 22d728fb7dSPeter Brune PetscFunctionBegin; 23d728fb7dSPeter Brune ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr); 24d728fb7dSPeter Brune snes = (SNES)ctx; 25be95d8f1SBarry Smith ierr = SNESGetNPC(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) { 29302440fdSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 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); 44724f369dSPeter Brune ierr = SNESGetJacobian(subsnes[i],&subJ,&subpJ,NULL,NULL);CHKERRQ(ierr); 45d728fb7dSPeter Brune ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr); 46724f369dSPeter Brune ierr = KSPSetOperators(ksp,subJ,subpJ);CHKERRQ(ierr); 47724f369dSPeter Brune ierr = KSPSolve(ksp,b[i],x[i]);CHKERRQ(ierr); 48d728fb7dSPeter Brune ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 49ed3c11a9SPeter Brune } 50ed3c11a9SPeter Brune for (i=0;i<n;i++) { 51d728fb7dSPeter Brune ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 52d728fb7dSPeter Brune } 53d728fb7dSPeter Brune PetscFunctionReturn(0); 54d728fb7dSPeter Brune } 55d728fb7dSPeter Brune 56d728fb7dSPeter Brune #undef __FUNCT__ 576cf69a2aSStefano Zampini #define __FUNCT__ "SNESDestroy_ASPIN" 586cf69a2aSStefano Zampini PetscErrorCode SNESDestroy_ASPIN(SNES snes) 596cf69a2aSStefano Zampini { 606cf69a2aSStefano Zampini PetscErrorCode ierr; 616cf69a2aSStefano Zampini 626cf69a2aSStefano Zampini PetscFunctionBegin; 636cf69a2aSStefano Zampini ierr = SNESDestroy(&snes->pc);CHKERRQ(ierr); 64263b3d60SStefano Zampini /* reset NEWTONLS and free the data */ 65263b3d60SStefano Zampini ierr = SNESReset(snes);CHKERRQ(ierr); 66263b3d60SStefano Zampini ierr = PetscFree(snes->data);CHKERRQ(ierr); 676cf69a2aSStefano Zampini PetscFunctionReturn(0); 686cf69a2aSStefano Zampini } 696cf69a2aSStefano Zampini 706cf69a2aSStefano Zampini #undef __FUNCT__ 71d728fb7dSPeter Brune #define __FUNCT__ "SNESCreate_ASPIN" 722e483eabSPeter Brune /* -------------------------------------------------------------------------- */ 732e483eabSPeter Brune /*MC 742e483eabSPeter Brune SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton 752e483eabSPeter Brune 762e483eabSPeter Brune Options Database: 772e483eabSPeter Brune + -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM) 782e483eabSPeter Brune . -npc_sub_snes_ - options prefix of the subdomain nonlinear solves 792e483eabSPeter Brune . -npc_sub_ksp_ - options prefix of the subdomain Krylov solver 802e483eabSPeter Brune - -npc_sub_pc_ - options prefix of the subdomain preconditioner 812e483eabSPeter Brune 822e483eabSPeter Brune Notes: This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other 832e483eabSPeter Brune similar functionality in SNES as it creates a linear shell matrix that corresponds to the product: 842e483eabSPeter Brune 852e483eabSPeter Brune \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b)) 862e483eabSPeter Brune 872e483eabSPeter Brune which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of 882e483eabSPeter Brune nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner 892e483eabSPeter Brune factorizations are reused on each application of J_b^{-1}. 902e483eabSPeter Brune 912e483eabSPeter Brune Level: intermediate 922e483eabSPeter Brune 93*4f02bc6aSBarry Smith References: 94*4f02bc6aSBarry Smith 95*4f02bc6aSBarry Smith X.-C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact {Newton} algorithms", SIAM J. Sci. Comput., 24, 2002. 96*4f02bc6aSBarry Smith 97*4f02bc6aSBarry Smith Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 98*4f02bc6aSBarry Smith SIAM Review, 57(4), 2015 99*4f02bc6aSBarry Smith 100be95d8f1SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetNPC(), SNESGetNPCSide() 1012e483eabSPeter Brune 1022e483eabSPeter Brune M*/ 1038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes) 104d728fb7dSPeter Brune { 105d728fb7dSPeter Brune PetscErrorCode ierr; 106d728fb7dSPeter Brune SNES npc; 107d728fb7dSPeter Brune KSP ksp; 108d728fb7dSPeter Brune PC pc; 109d728fb7dSPeter Brune Mat aspinmat; 110d728fb7dSPeter Brune Vec F; 111d728fb7dSPeter Brune PetscInt n; 112ed3c11a9SPeter Brune SNESLineSearch linesearch; 113d728fb7dSPeter Brune 114d728fb7dSPeter Brune PetscFunctionBegin; 115d728fb7dSPeter Brune /* set up the solver */ 116d728fb7dSPeter Brune ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr); 117be95d8f1SBarry Smith ierr = SNESSetNPCSide(snes,PC_LEFT);CHKERRQ(ierr); 118ce5a860aSPeter Brune ierr = SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED);CHKERRQ(ierr); 119be95d8f1SBarry Smith ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr); 120d728fb7dSPeter Brune ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr); 121d78e9beeSPeter Brune ierr = SNESNASMSetType(npc,PC_ASM_BASIC);CHKERRQ(ierr); 122d728fb7dSPeter Brune ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr); 123d728fb7dSPeter Brune ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 124d728fb7dSPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 125d728fb7dSPeter Brune ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1267601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 127ed3c11a9SPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);CHKERRQ(ierr); 128d728fb7dSPeter Brune 129d728fb7dSPeter Brune /* set up the shell matrix */ 130d728fb7dSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 131d728fb7dSPeter Brune ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr); 1326cf69a2aSStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)snes),n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr); 133d728fb7dSPeter Brune ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr); 134d728fb7dSPeter Brune ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr); 135d728fb7dSPeter Brune ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr); 136d728fb7dSPeter Brune ierr = MatDestroy(&aspinmat);CHKERRQ(ierr); 137d728fb7dSPeter Brune 1386cf69a2aSStefano Zampini snes->ops->destroy = SNESDestroy_ASPIN; 1396cf69a2aSStefano Zampini 140d728fb7dSPeter Brune PetscFunctionReturn(0); 141d728fb7dSPeter Brune } 142