Actual source code: snesmfj.c
petsc-3.3-p7 2013-05-11
2: #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
3: #include <../src/mat/impls/mffd/mffdimpl.h>
4: #include <petsc-private/matimpl.h>
8: /*@C
9: MatMFFDComputeJacobian - Tells the matrix-free Jacobian object the new location at which
10: Jacobian matrix vector products will be computed at, i.e. J(x) * a. The x is obtained
11: from the SNES object (using SNESGetSolution()).
13: Logically Collective on SNES
15: Input Parameters:
16: + snes - the nonlinear solver context
17: . x - the point at which the Jacobian vector products will be performed
18: . jac - the matrix-free Jacobian object
19: . B - either the same as jac or another matrix type (ignored)
20: . flag - not relevent for matrix-free form
21: - dummy - the user context (ignored)
23: Level: developer
25: Warning:
26: If MatMFFDSetBase() is ever called on jac then this routine will NO longer get
27: the x from the SNES object and MatMFFDSetBase() must from that point on be used to
28: change the base vector x.
30: Notes:
31: This can be passed into SNESSetJacobian() as the Jacobian evaluation function argument
32: when using a completely matrix-free solver,
33: that is the B matrix is also the same matrix operator. This is used when you select
34: -snes_mf but rarely used directly by users. (All this routine does is call MatAssemblyBegin/End() on
35: the Mat jac.
37: .seealso: MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD,
38: MatMFFDSetHHistory(), MatMFFDSetFunctionError(), MatCreateMFFD(), SNESSetJacobian()
40: @*/
41: PetscErrorCode MatMFFDComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
42: {
45: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
47: return(0);
48: }
50: PetscErrorCode MatAssemblyEnd_MFFD(Mat,MatAssemblyType);
51: EXTERN_C_BEGIN
52: PetscErrorCode MatMFFDSetBase_MFFD(Mat,Vec,Vec);
53: EXTERN_C_END
57: /*
58: MatAssemblyEnd_SNESMF - Calls MatAssemblyEnd_MFFD() and then sets the
59: base from the SNES context
61: */
62: PetscErrorCode MatAssemblyEnd_SNESMF(Mat J,MatAssemblyType mt)
63: {
65: MatMFFD j = (MatMFFD)J->data;
66: SNES snes = (SNES)j->funcctx;
67: Vec u,f;
70: MatAssemblyEnd_MFFD(J,mt);
72: SNESGetSolution(snes,&u);
73: SNESGetFunction(snes,&f,PETSC_NULL,PETSC_NULL);
74: MatMFFDSetBase_MFFD(J,u,f);
75: return(0);
76: }
78: EXTERN_C_BEGIN
79: /*
80: This routine resets the MatAssemblyEnd() for the MatMFFD created from MatCreateSNESMF() so that it NO longer
81: uses the solution in the SNES object to update the base. See the warning in MatCreateSNESMF().
82: */
85: PetscErrorCode MatMFFDSetBase_SNESMF(Mat J,Vec U,Vec F)
86: {
90: MatMFFDSetBase_MFFD(J,U,F);
91: J->ops->assemblyend = MatAssemblyEnd_MFFD;
92: return(0);
93: }
94: EXTERN_C_END
98: /*@
99: MatCreateSNESMF - Creates a matrix-free matrix context for use with
100: a SNES solver. This matrix can be used as the Jacobian argument for
101: the routine SNESSetJacobian(). See MatCreateMFFD() for details on how
102: the finite difference computation is done.
104: Collective on SNES and Vec
106: Input Parameters:
107: . snes - the SNES context
109: Output Parameter:
110: . J - the matrix-free matrix
112: Level: advanced
114: Warning:
115: If MatMFFDSetBase() is ever called on jac then this routine will NO longer get
116: the x from the SNES object and MatMFFDSetBase() must from that point on be used to
117: change the base vector x.
119: Notes: The difference between this routine and MatCreateMFFD() is that this matrix
120: automatically gets the current base vector from the SNES object and not from an
121: explicit call to MatMFFDSetBase().
123: .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin()
124: MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateMFFD(),
125: MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian()
126:
127: @*/
128: PetscErrorCode MatCreateSNESMF(SNES snes,Mat *J)
129: {
131: PetscInt n,N;
134: if (snes->vec_func) {
135: VecGetLocalSize(snes->vec_func,&n);
136: VecGetSize(snes->vec_func,&N);
137: } else if (snes->dm) {
138: Vec tmp;
139: DMGetGlobalVector(snes->dm,&tmp);
140: VecGetLocalSize(tmp,&n);
141: VecGetSize(tmp,&N);
142: DMRestoreGlobalVector(snes->dm,&tmp);
143: } else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
144: MatCreateMFFD(((PetscObject)snes)->comm,n,n,N,N,J);
145: MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
146: (*J)->ops->assemblyend = MatAssemblyEnd_SNESMF;
147: PetscObjectComposeFunctionDynamic((PetscObject)*J,"MatMFFDSetBase_C","MatMFFDSetBase_SNESMF",MatMFFDSetBase_SNESMF);
148: return(0);
149: }