1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2d728fb7dSPeter Brune #include <petscdm.h> 3d728fb7dSPeter Brune 4d728fb7dSPeter Brune PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y) 5d728fb7dSPeter Brune { 6d728fb7dSPeter Brune PetscErrorCode ierr; 7d728fb7dSPeter Brune void *ctx; 8d728fb7dSPeter Brune SNES snes; 9d728fb7dSPeter Brune PetscInt n,i; 10d728fb7dSPeter Brune VecScatter *oscatter; 11d728fb7dSPeter Brune SNES *subsnes; 12d728fb7dSPeter Brune PetscBool match; 13d728fb7dSPeter Brune MPI_Comm comm; 14d728fb7dSPeter Brune KSP ksp; 15d728fb7dSPeter Brune Vec *x,*b; 16ed3c11a9SPeter Brune Vec W; 17d728fb7dSPeter Brune SNES npc; 18724f369dSPeter Brune Mat subJ,subpJ; 19d728fb7dSPeter Brune 20d728fb7dSPeter Brune PetscFunctionBegin; 21d728fb7dSPeter Brune ierr = MatShellGetContext(m,&ctx);CHKERRQ(ierr); 22d728fb7dSPeter Brune snes = (SNES)ctx; 23be95d8f1SBarry Smith ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr); 24ed3c11a9SPeter Brune ierr = SNESGetFunction(npc,&W,NULL,NULL);CHKERRQ(ierr); 25d728fb7dSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);CHKERRQ(ierr); 26d728fb7dSPeter Brune if (!match) { 27302440fdSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 28d728fb7dSPeter Brune SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz"); 29d728fb7dSPeter Brune } 30d728fb7dSPeter Brune ierr = SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);CHKERRQ(ierr); 31d728fb7dSPeter Brune ierr = SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);CHKERRQ(ierr); 32ed3c11a9SPeter Brune 33ed3c11a9SPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 34ed3c11a9SPeter Brune ierr = MatMult(npc->jacobian_pre,X,W);CHKERRQ(ierr); 35ed3c11a9SPeter Brune 36d728fb7dSPeter Brune for (i=0;i<n;i++) { 37ed3c11a9SPeter Brune ierr = VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 38d728fb7dSPeter Brune } 39d728fb7dSPeter Brune for (i=0;i<n;i++) { 40ed3c11a9SPeter Brune ierr = VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 41ed3c11a9SPeter Brune ierr = VecSet(x[i],0.);CHKERRQ(ierr); 42724f369dSPeter Brune ierr = SNESGetJacobian(subsnes[i],&subJ,&subpJ,NULL,NULL);CHKERRQ(ierr); 43d728fb7dSPeter Brune ierr = SNESGetKSP(subsnes[i],&ksp);CHKERRQ(ierr); 44724f369dSPeter Brune ierr = KSPSetOperators(ksp,subJ,subpJ);CHKERRQ(ierr); 45724f369dSPeter Brune ierr = KSPSolve(ksp,b[i],x[i]);CHKERRQ(ierr); 46d728fb7dSPeter Brune ierr = VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 47d728fb7dSPeter Brune ierr = VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 48d728fb7dSPeter Brune } 49d728fb7dSPeter Brune PetscFunctionReturn(0); 50d728fb7dSPeter Brune } 51d728fb7dSPeter Brune 52e0877f53SBarry Smith static PetscErrorCode SNESDestroy_ASPIN(SNES snes) 536cf69a2aSStefano Zampini { 546cf69a2aSStefano Zampini PetscErrorCode ierr; 556cf69a2aSStefano Zampini 566cf69a2aSStefano Zampini PetscFunctionBegin; 57efd4aadfSBarry Smith ierr = SNESDestroy(&snes->npc);CHKERRQ(ierr); 58263b3d60SStefano Zampini /* reset NEWTONLS and free the data */ 59263b3d60SStefano Zampini ierr = SNESReset(snes);CHKERRQ(ierr); 60263b3d60SStefano Zampini ierr = PetscFree(snes->data);CHKERRQ(ierr); 616cf69a2aSStefano Zampini PetscFunctionReturn(0); 626cf69a2aSStefano Zampini } 636cf69a2aSStefano Zampini 642e483eabSPeter Brune /* -------------------------------------------------------------------------- */ 652e483eabSPeter Brune /*MC 662e483eabSPeter Brune SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton 672e483eabSPeter Brune 682e483eabSPeter Brune Options Database: 692e483eabSPeter Brune + -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM) 702e483eabSPeter Brune . -npc_sub_snes_ - options prefix of the subdomain nonlinear solves 712e483eabSPeter Brune . -npc_sub_ksp_ - options prefix of the subdomain Krylov solver 722e483eabSPeter Brune - -npc_sub_pc_ - options prefix of the subdomain preconditioner 732e483eabSPeter Brune 74d14a3f37SRichard Tran Mills Notes: 75d14a3f37SRichard Tran Mills This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other 76d14a3f37SRichard Tran Mills similar functionality in SNES as it creates a linear shell matrix that corresponds to the product 772e483eabSPeter Brune 782e483eabSPeter Brune \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b)) 792e483eabSPeter Brune 802e483eabSPeter Brune which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of 812e483eabSPeter Brune nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner 822e483eabSPeter Brune factorizations are reused on each application of J_b^{-1}. 832e483eabSPeter Brune 84951fe5abSBarry Smith The Krylov method used in this nonlinear solver is run with NO preconditioner, because the preconditioning is done 85951fe5abSBarry Smith at the nonlinear level, but the Jacobian for the original function must be provided (or calculated via coloring and 86951fe5abSBarry Smith finite differences automatically) in the Pmat location of SNESSetJacobian() because the action of the original Jacobian 87951fe5abSBarry Smith is needed by the shell matrix used to apply the Jacobian of the nonlinear preconditioned problem (see above). 88951fe5abSBarry Smith Note that since the Pmat is not used to construct a preconditioner it could be provided in a matrix-free form. 89951fe5abSBarry Smith The code for this implementation is a bit confusing because the Amat of SNESSetJacobian() applies the Jacobian of the 90951fe5abSBarry Smith nonlinearly preconditioned function Jacobian while the Pmat provides the Jacobian of the original user provided function. 91951fe5abSBarry Smith Note that the original SNES and nonlinear preconditioner preconditioner (see SNESGetNPC()), in this case NASM, share 92e0c02066SRichard Tran Mills the same Jacobian matrices. SNESNASM computes the needed Jacobian in SNESNASMComputeFinalJacobian_Private(). 93951fe5abSBarry Smith 9406dd6b0eSSatish Balay Level: intermediate 952e483eabSPeter Brune 964f02bc6aSBarry Smith References: 97*606c0280SSatish Balay + * - X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms", SIAM J. Sci. Comput., 24, 2002. 98*606c0280SSatish Balay - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 994f02bc6aSBarry Smith SIAM Review, 57(4), 2015 1004f02bc6aSBarry Smith 101be95d8f1SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetNPC(), SNESGetNPCSide() 1022e483eabSPeter Brune 1032e483eabSPeter Brune M*/ 1048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes) 105d728fb7dSPeter Brune { 106d728fb7dSPeter Brune PetscErrorCode ierr; 107d728fb7dSPeter Brune SNES npc; 108d728fb7dSPeter Brune KSP ksp; 109d728fb7dSPeter Brune PC pc; 110d728fb7dSPeter Brune Mat aspinmat; 111d728fb7dSPeter Brune Vec F; 112d728fb7dSPeter Brune PetscInt n; 113ed3c11a9SPeter Brune SNESLineSearch linesearch; 114d728fb7dSPeter Brune 115d728fb7dSPeter Brune PetscFunctionBegin; 116d728fb7dSPeter Brune /* set up the solver */ 117d728fb7dSPeter Brune ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr); 118be95d8f1SBarry Smith ierr = SNESSetNPCSide(snes,PC_LEFT);CHKERRQ(ierr); 119ce5a860aSPeter Brune ierr = SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED);CHKERRQ(ierr); 120be95d8f1SBarry Smith ierr = SNESGetNPC(snes,&npc);CHKERRQ(ierr); 121d728fb7dSPeter Brune ierr = SNESSetType(npc,SNESNASM);CHKERRQ(ierr); 122d78e9beeSPeter Brune ierr = SNESNASMSetType(npc,PC_ASM_BASIC);CHKERRQ(ierr); 123d728fb7dSPeter Brune ierr = SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);CHKERRQ(ierr); 124d728fb7dSPeter Brune ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 125d728fb7dSPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 126d728fb7dSPeter Brune ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1277601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 128ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 129ed3c11a9SPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);CHKERRQ(ierr); 130ec786807SJed Brown } 131d728fb7dSPeter Brune 132d728fb7dSPeter Brune /* set up the shell matrix */ 133d728fb7dSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 134d728fb7dSPeter Brune ierr = VecGetLocalSize(F,&n);CHKERRQ(ierr); 1356cf69a2aSStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)snes),n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);CHKERRQ(ierr); 136d728fb7dSPeter Brune ierr = MatSetType(aspinmat,MATSHELL);CHKERRQ(ierr); 137d728fb7dSPeter Brune ierr = MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);CHKERRQ(ierr); 138d728fb7dSPeter Brune ierr = SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);CHKERRQ(ierr); 139d728fb7dSPeter Brune ierr = MatDestroy(&aspinmat);CHKERRQ(ierr); 140d728fb7dSPeter Brune 1416cf69a2aSStefano Zampini snes->ops->destroy = SNESDestroy_ASPIN; 1426cf69a2aSStefano Zampini 143d728fb7dSPeter Brune PetscFunctionReturn(0); 144d728fb7dSPeter Brune } 145