Actual source code: aspin.c
petsc-3.9.4 2018-09-11
1: #include <petsc/private/snesimpl.h>
2: #include <petscdm.h>
4: PetscErrorCode MatMultASPIN(Mat m,Vec X,Vec Y)
5: {
7: void *ctx;
8: SNES snes;
9: PetscInt n,i;
10: VecScatter *oscatter;
11: SNES *subsnes;
12: PetscBool match;
13: MPI_Comm comm;
14: KSP ksp;
15: Vec *x,*b;
16: Vec W;
17: SNES npc;
18: Mat subJ,subpJ;
21: MatShellGetContext(m,&ctx);
22: snes = (SNES)ctx;
23: SNESGetNPC(snes,&npc);
24: SNESGetFunction(npc,&W,NULL,NULL);
25: PetscObjectTypeCompare((PetscObject)npc,SNESNASM,&match);
26: if (!match) {
27: PetscObjectGetComm((PetscObject)snes,&comm);
28: SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"MatMultASPIN requires that the nonlinear preconditioner be Nonlinear additive Schwarz");
29: }
30: SNESNASMGetSubdomains(npc,&n,&subsnes,NULL,&oscatter,NULL);
31: SNESNASMGetSubdomainVecs(npc,&n,&x,&b,NULL,NULL);
33: VecSet(Y,0);
34: MatMult(npc->jacobian_pre,X,W);
36: for (i=0;i<n;i++) {
37: VecScatterBegin(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);
38: }
39: for (i=0;i<n;i++) {
40: VecScatterEnd(oscatter[i],W,b[i],INSERT_VALUES,SCATTER_FORWARD);
41: VecSet(x[i],0.);
42: SNESGetJacobian(subsnes[i],&subJ,&subpJ,NULL,NULL);
43: SNESGetKSP(subsnes[i],&ksp);
44: KSPSetOperators(ksp,subJ,subpJ);
45: KSPSolve(ksp,b[i],x[i]);
46: VecScatterBegin(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);
47: }
48: for (i=0;i<n;i++) {
49: VecScatterEnd(oscatter[i],x[i],Y,ADD_VALUES,SCATTER_REVERSE);
50: }
51: return(0);
52: }
54: static PetscErrorCode SNESDestroy_ASPIN(SNES snes)
55: {
59: SNESDestroy(&snes->npc);
60: /* reset NEWTONLS and free the data */
61: SNESReset(snes);
62: PetscFree(snes->data);
63: return(0);
64: }
66: /* -------------------------------------------------------------------------- */
67: /*MC
68: SNESASPIN - Helper SNES type for Additive-Schwarz Preconditioned Inexact Newton
70: Options Database:
71: + -npc_snes_ - options prefix of the nonlinear subdomain solver (must be of type NASM)
72: . -npc_sub_snes_ - options prefix of the subdomain nonlinear solves
73: . -npc_sub_ksp_ - options prefix of the subdomain Krylov solver
74: - -npc_sub_pc_ - options prefix of the subdomain preconditioner
76: Notes:
77: This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other
78: similar functionality in SNES as it creates a linear shell matrix that corresponds to the product
80: \sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
82: which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
83: nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
84: factorizations are reused on each application of J_b^{-1}.
86: Level:
87: intermediate
89: References:
90: + 1. - X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms", SIAM J. Sci. Comput., 24, 2002.
91: - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
92: SIAM Review, 57(4), 2015
94: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNASM, SNESGetNPC(), SNESGetNPCSide()
96: M*/
97: PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES snes)
98: {
100: SNES npc;
101: KSP ksp;
102: PC pc;
103: Mat aspinmat;
104: Vec F;
105: PetscInt n;
106: SNESLineSearch linesearch;
109: /* set up the solver */
110: SNESSetType(snes,SNESNEWTONLS);
111: SNESSetNPCSide(snes,PC_LEFT);
112: SNESSetFunctionType(snes,SNES_FUNCTION_PRECONDITIONED);
113: SNESGetNPC(snes,&npc);
114: SNESSetType(npc,SNESNASM);
115: SNESNASMSetType(npc,PC_ASM_BASIC);
116: SNESNASMSetComputeFinalJacobian(npc,PETSC_TRUE);
117: SNESGetKSP(snes,&ksp);
118: KSPGetPC(ksp,&pc);
119: PCSetType(pc,PCNONE);
120: SNESGetLineSearch(snes,&linesearch);
121: SNESLineSearchSetType(linesearch,SNESLINESEARCHBT);
123: /* set up the shell matrix */
124: SNESGetFunction(snes,&F,NULL,NULL);
125: VecGetLocalSize(F,&n);
126: MatCreateShell(PetscObjectComm((PetscObject)snes),n,n,PETSC_DECIDE,PETSC_DECIDE,snes,&aspinmat);
127: MatSetType(aspinmat,MATSHELL);
128: MatShellSetOperation(aspinmat,MATOP_MULT,(void(*)(void))MatMultASPIN);
129: SNESSetJacobian(snes,aspinmat,NULL,NULL,NULL);
130: MatDestroy(&aspinmat);
132: snes->ops->destroy = SNESDestroy_ASPIN;
134: return(0);
135: }