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 void *ctx; 7d728fb7dSPeter Brune SNES snes; 8d728fb7dSPeter Brune PetscInt n,i; 9d728fb7dSPeter Brune VecScatter *oscatter; 10d728fb7dSPeter Brune SNES *subsnes; 11d728fb7dSPeter Brune PetscBool match; 12d728fb7dSPeter Brune MPI_Comm comm; 13d728fb7dSPeter Brune KSP ksp; 14d728fb7dSPeter Brune Vec *x,*b; 15ed3c11a9SPeter Brune Vec W; 16d728fb7dSPeter Brune SNES npc; 17724f369dSPeter Brune Mat subJ,subpJ; 18d728fb7dSPeter Brune 19d728fb7dSPeter Brune PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(m,&ctx)); 21d728fb7dSPeter Brune snes = (SNES)ctx; 229566063dSJacob Faibussowitsch PetscCall(SNESGetNPC(snes,&npc)); 239566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(npc,&W,NULL,NULL)); 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match)); 25d728fb7dSPeter Brune if (!match) { 269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 27d728fb7dSPeter Brune SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz"); 28d728fb7dSPeter Brune } 299566063dSJacob Faibussowitsch PetscCall(SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL)); 309566063dSJacob Faibussowitsch PetscCall(SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL)); 31ed3c11a9SPeter Brune 329566063dSJacob Faibussowitsch PetscCall(VecSet(Y,0)); 339566063dSJacob Faibussowitsch PetscCall(MatMult(npc->jacobian_pre,X,W)); 34ed3c11a9SPeter Brune 35d728fb7dSPeter Brune for (i=0;i<n;i++) { 369566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD)); 37d728fb7dSPeter Brune } 38d728fb7dSPeter Brune for (i=0;i<n;i++) { 399566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD)); 409566063dSJacob Faibussowitsch PetscCall(VecSet(x[i],0.)); 419566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(subsnes[i],&subJ,&subpJ,NULL,NULL)); 429566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(subsnes[i],&ksp)); 439566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp,subJ,subpJ)); 449566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp,b[i],x[i])); 459566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE)); 469566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE)); 47d728fb7dSPeter Brune } 48d728fb7dSPeter Brune PetscFunctionReturn(0); 49d728fb7dSPeter Brune } 50d728fb7dSPeter Brune 51e0877f53SBarry Smith static PetscErrorCode SNESDestroy_ASPIN(SNES snes) 526cf69a2aSStefano Zampini { 536cf69a2aSStefano Zampini PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes->npc)); 55263b3d60SStefano Zampini /* reset NEWTONLS and free the data */ 569566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 579566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 586cf69a2aSStefano Zampini PetscFunctionReturn(0); 596cf69a2aSStefano Zampini } 606cf69a2aSStefano Zampini 612e483eabSPeter Brune /* -------------------------------------------------------------------------- */ 622e483eabSPeter Brune /*MC 632e483eabSPeter Brune SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton 642e483eabSPeter Brune 652e483eabSPeter Brune Options Database: 662e483eabSPeter Brune + -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM) 672e483eabSPeter Brune . -npc_sub_snes_ - options prefix of the subdomain nonlinear solves 682e483eabSPeter Brune . -npc_sub_ksp_ - options prefix of the subdomain Krylov solver 692e483eabSPeter Brune - -npc_sub_pc_ - options prefix of the subdomain preconditioner 702e483eabSPeter Brune 71d14a3f37SRichard Tran Mills Notes: 72d14a3f37SRichard Tran Mills This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other 73d14a3f37SRichard Tran Mills similar functionality in SNES as it creates a linear shell matrix that corresponds to the product 742e483eabSPeter Brune 752e483eabSPeter Brune \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b)) 762e483eabSPeter Brune 772e483eabSPeter Brune which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of 782e483eabSPeter Brune nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner 792e483eabSPeter Brune factorizations are reused on each application of J_b^{-1}. 802e483eabSPeter Brune 81951fe5abSBarry Smith The Krylov method used in this nonlinear solver is run with NO preconditioner, because the preconditioning is done 82951fe5abSBarry Smith at the nonlinear level, but the Jacobian for the original function must be provided (or calculated via coloring and 83951fe5abSBarry Smith finite differences automatically) in the Pmat location of SNESSetJacobian() because the action of the original Jacobian 84951fe5abSBarry Smith is needed by the shell matrix used to apply the Jacobian of the nonlinear preconditioned problem (see above). 85951fe5abSBarry Smith Note that since the Pmat is not used to construct a preconditioner it could be provided in a matrix-free form. 86951fe5abSBarry Smith The code for this implementation is a bit confusing because the Amat of SNESSetJacobian() applies the Jacobian of the 87951fe5abSBarry Smith nonlinearly preconditioned function Jacobian while the Pmat provides the Jacobian of the original user provided function. 88951fe5abSBarry Smith Note that the original SNES and nonlinear preconditioner preconditioner (see SNESGetNPC()), in this case NASM, share 89e0c02066SRichard Tran Mills the same Jacobian matrices. SNESNASM computes the needed Jacobian in SNESNASMComputeFinalJacobian_Private(). 90951fe5abSBarry Smith 9106dd6b0eSSatish Balay Level: intermediate 922e483eabSPeter Brune 934f02bc6aSBarry Smith References: 94606c0280SSatish Balay + * - X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms", SIAM J. Sci. Comput., 24, 2002. 95606c0280SSatish Balay - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 964f02bc6aSBarry Smith SIAM Review, 57(4), 2015 974f02bc6aSBarry Smith 98*db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNASM`, `SNESGetNPC()`, `SNESGetNPCSide()` 992e483eabSPeter Brune 1002e483eabSPeter Brune M*/ 1018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes) 102d728fb7dSPeter Brune { 103d728fb7dSPeter Brune SNES npc; 104d728fb7dSPeter Brune KSP ksp; 105d728fb7dSPeter Brune PC pc; 106d728fb7dSPeter Brune Mat aspinmat; 107d728fb7dSPeter Brune Vec F; 108d728fb7dSPeter Brune PetscInt n; 109ed3c11a9SPeter Brune SNESLineSearch linesearch; 110d728fb7dSPeter Brune 111d728fb7dSPeter Brune PetscFunctionBegin; 112d728fb7dSPeter Brune /* set up the solver */ 1139566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes,SNESNEWTONLS)); 1149566063dSJacob Faibussowitsch PetscCall(SNESSetNPCSide(snes,PC_LEFT)); 1159566063dSJacob Faibussowitsch PetscCall(SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED)); 1169566063dSJacob Faibussowitsch PetscCall(SNESGetNPC(snes,&npc)); 1179566063dSJacob Faibussowitsch PetscCall(SNESSetType(npc,SNESNASM)); 1189566063dSJacob Faibussowitsch PetscCall(SNESNASMSetType(npc,PC_ASM_BASIC)); 1199566063dSJacob Faibussowitsch PetscCall(SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE)); 1209566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 1219566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 1229566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 1239566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes,&linesearch)); 124ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 1259566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBT)); 126ec786807SJed Brown } 127d728fb7dSPeter Brune 128d728fb7dSPeter Brune /* set up the shell matrix */ 1299566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,&F,NULL,NULL)); 1309566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(F,&n)); 1319566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)snes),n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetType(aspinmat,MATSHELL)); 1339566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN)); 1349566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL)); 1359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aspinmat)); 136d728fb7dSPeter Brune 1376cf69a2aSStefano Zampini snes->ops->destroy = SNESDestroy_ASPIN; 1386cf69a2aSStefano Zampini 139d728fb7dSPeter Brune PetscFunctionReturn(0); 140d728fb7dSPeter Brune } 141