Actual source code: snesmfj.c
petsc-3.4.2 2013-07-02
2: #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
3: #include <petscdm.h> /*I "petscdm.h" I*/
4: #include <../src/mat/impls/mffd/mffdimpl.h>
5: #include <petsc-private/matimpl.h>
9: /*@C
10: MatMFFDComputeJacobian - Tells the matrix-free Jacobian object the new location at which
11: Jacobian matrix vector products will be computed at, i.e. J(x) * a. The x is obtained
12: from the SNES object (using SNESGetSolution()).
14: Logically Collective on SNES
16: Input Parameters:
17: + snes - the nonlinear solver context
18: . x - the point at which the Jacobian vector products will be performed
19: . jac - the matrix-free Jacobian object
20: . B - either the same as jac or another matrix type (ignored)
21: . flag - not relevent for matrix-free form
22: - dummy - the user context (ignored)
24: Level: developer
26: Warning:
27: If MatMFFDSetBase() is ever called on jac then this routine will NO longer get
28: the x from the SNES object and MatMFFDSetBase() must from that point on be used to
29: change the base vector x.
31: Notes:
32: This can be passed into SNESSetJacobian() as the Jacobian evaluation function argument
33: when using a completely matrix-free solver,
34: that is the B matrix is also the same matrix operator. This is used when you select
35: -snes_mf but rarely used directly by users. (All this routine does is call MatAssemblyBegin/End() on
36: the Mat jac.
38: .seealso: MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD,
39: MatMFFDSetHHistory(), MatMFFDSetFunctionError(), MatCreateMFFD(), SNESSetJacobian()
41: @*/
42: PetscErrorCode MatMFFDComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
43: {
47: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
48: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
49: return(0);
50: }
52: PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat,MatAssemblyType);
53: PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat,Vec,Vec);
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,NULL,NULL);
74: MatMFFDSetBase_MFFD(J,u,f);
75: return(0);
76: }
78: /*
79: This routine resets the MatAssemblyEnd() for the MatMFFD created from MatCreateSNESMF() so that it NO longer
80: uses the solution in the SNES object to update the base. See the warning in MatCreateSNESMF().
81: */
84: PetscErrorCode MatMFFDSetBase_SNESMF(Mat J,Vec U,Vec F)
85: {
89: MatMFFDSetBase_MFFD(J,U,F);
91: J->ops->assemblyend = MatAssemblyEnd_MFFD;
92: return(0);
93: }
97: /*@
98: MatCreateSNESMF - Creates a matrix-free matrix context for use with
99: a SNES solver. This matrix can be used as the Jacobian argument for
100: the routine SNESSetJacobian(). See MatCreateMFFD() for details on how
101: the finite difference computation is done.
103: Collective on SNES and Vec
105: Input Parameters:
106: . snes - the SNES context
108: Output Parameter:
109: . J - the matrix-free matrix
111: Level: advanced
113: Warning:
114: If MatMFFDSetBase() is ever called on jac then this routine will NO longer get
115: the x from the SNES object and MatMFFDSetBase() must from that point on be used to
116: change the base vector x.
118: Notes: The difference between this routine and MatCreateMFFD() is that this matrix
119: automatically gets the current base vector from the SNES object and not from an
120: explicit call to MatMFFDSetBase().
122: .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin()
123: MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateMFFD(),
124: MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
126: @*/
127: PetscErrorCode MatCreateSNESMF(SNES snes,Mat *J)
128: {
130: PetscInt n,N;
133: if (snes->vec_func) {
134: VecGetLocalSize(snes->vec_func,&n);
135: VecGetSize(snes->vec_func,&N);
136: } else if (snes->dm) {
137: Vec tmp;
138: DMGetGlobalVector(snes->dm,&tmp);
139: VecGetLocalSize(tmp,&n);
140: VecGetSize(tmp,&N);
141: DMRestoreGlobalVector(snes->dm,&tmp);
142: } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
143: MatCreateMFFD(PetscObjectComm((PetscObject)snes),n,n,N,N,J);
144: MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
146: (*J)->ops->assemblyend = MatAssemblyEnd_SNESMF;
148: PetscObjectComposeFunction((PetscObject)*J,"MatMFFDSetBase_C",MatMFFDSetBase_SNESMF);
149: return(0);
150: }