Actual source code: snesmfj.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <petsc/private/snesimpl.h>
  3:  #include <petscdm.h>
  4:  #include <../src/mat/impls/mffd/mffdimpl.h>
  5:  #include <petsc/private/matimpl.h>

  7: /*@C
  8:    MatMFFDComputeJacobian - Tells the matrix-free Jacobian object the new location at which
  9:        Jacobian matrix vector products will be computed at, i.e. J(x) * a. The x is obtained
 10:        from the SNES object (using SNESGetSolution()).

 12:    Logically Collective on SNES

 14:    Input Parameters:
 15: +   snes - the nonlinear solver context
 16: .   x - the point at which the Jacobian vector products will be performed
 17: .   jac - the matrix-free Jacobian object
 18: .   B - either the same as jac or another matrix type (ignored)
 19: .   flag - not relevent for matrix-free form
 20: -   dummy - the user context (ignored)

 22:    Level: developer

 24:    Warning:
 25:       If MatMFFDSetBase() is ever called on jac then this routine will NO longer get
 26:     the x from the SNES object and MatMFFDSetBase() must from that point on be used to
 27:     change the base vector x.

 29:    Notes:
 30:      This can be passed into SNESSetJacobian() as the Jacobian evaluation function argument
 31:      when using a completely matrix-free solver,
 32:      that is the B matrix is also the same matrix operator. This is used when you select
 33:      -snes_mf but rarely used directly by users. (All this routine does is call MatAssemblyBegin/End() on
 34:      the Mat jac.)

 36: .seealso: MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD,
 37:           MatMFFDSetHHistory(), MatMFFDSetFunctionError(), MatCreateMFFD(), SNESSetJacobian()

 39: @*/
 40: PetscErrorCode  MatMFFDComputeJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
 41: {

 45:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
 46:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
 47:   return(0);
 48: }

 50: PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat,MatAssemblyType);
 51: PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat,Vec,Vec);

 53: /*
 54:    MatAssemblyEnd_SNESMF - Calls MatAssemblyEnd_MFFD() and then sets the
 55:     base from the SNES context

 57: */
 58: static PetscErrorCode MatAssemblyEnd_SNESMF(Mat J,MatAssemblyType mt)
 59: {
 61:   MatMFFD        j    = (MatMFFD)J->data;
 62:   SNES           snes = (SNES)j->ctx;
 63:   Vec            u,f;

 66:   MatAssemblyEnd_MFFD(J,mt);

 68:   SNESGetSolution(snes,&u);
 69:   if (j->func == (PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction) {
 70:     SNESGetFunction(snes,&f,NULL,NULL);
 71:     MatMFFDSetBase_MFFD(J,u,f);
 72:   } else {
 73:     /* f value known by SNES is not correct for other differencing function */
 74:     MatMFFDSetBase_MFFD(J,u,NULL);
 75:   }
 76:   return(0);
 77: }

 79: /*
 80:    MatAssemblyEnd_SNESMF_UseBase - Calls MatAssemblyEnd_MFFD() and then sets the
 81:     base from the SNES context. This version will cause the base to be used for differencing
 82:     even if the func is not SNESComputeFunction. See: MatSNESMFUseBase()


 85: */
 86: static PetscErrorCode MatAssemblyEnd_SNESMF_UseBase(Mat J,MatAssemblyType mt)
 87: {
 89:   MatMFFD        j    = (MatMFFD)J->data;
 90:   SNES           snes = (SNES)j->ctx;
 91:   Vec            u,f;

 94:   MatAssemblyEnd_MFFD(J,mt);

 96:   SNESGetSolution(snes,&u);
 97:   SNESGetFunction(snes,&f,NULL,NULL);
 98:   MatMFFDSetBase_MFFD(J,u,f);
 99:   return(0);
100: }

102: /*
103:     This routine resets the MatAssemblyEnd() for the MatMFFD created from MatCreateSNESMF() so that it NO longer
104:   uses the solution in the SNES object to update the base. See the warning in MatCreateSNESMF().
105: */
106: static PetscErrorCode  MatMFFDSetBase_SNESMF(Mat J,Vec U,Vec F)
107: {

111:   MatMFFDSetBase_MFFD(J,U,F);
112:   J->ops->assemblyend = MatAssemblyEnd_MFFD;
113:   return(0);
114: }

116: static PetscErrorCode  MatSNESMFSetReuseBase_SNESMF(Mat J,PetscBool use)
117: {
119:   if (use) {
120:     J->ops->assemblyend = MatAssemblyEnd_SNESMF_UseBase;
121:   } else {
122:     J->ops->assemblyend = MatAssemblyEnd_SNESMF;
123:   }
124:   return(0);
125: }

127: /*@
128:     MatSNESMFSetReuseBase - Causes the base vector to be used for differencing even if the function provided to SNESSetFunction() is not the
129:                        same as that provided to MatMFFDSetFunction().

131:     Logically Collective on Mat

133:     Input Parameters:
134: +   J - the MatMFFD matrix
135: -   use - if true always reuse the base vector instead of recomputing f(u) even if the function in the MatSNESMF is 
136:           not SNESComputeFunction()

138:     Notes: Care must be taken when using this routine to insure that the function provided to MatMFFDSetFunction(), call it F_MF() is compatible with
139:            with that provided to SNESSetFunction(), call it F_SNES(). That is, (F_MF(u + h*d) - F_SNES(u))/h has to approximate the derivative

141:     Developer Notes: This was provided for the MOOSE team who desired to have a SNESSetFunction() function that could change configurations due
142:                      to contacts while the function provided to MatMFFDSetFunction() cannot. Except for the possibility of changing the configuration
143:                      both functions compute the same mathematical function so the differencing makes sense.

145:     Level: advanced

147: .seealso: MatCreateSNESMF(), MatSNESMFGetReuseBase()
148: @*/
149: PetscErrorCode  MatSNESMFSetReuseBase(Mat J,PetscBool use)
150: {

155:   PetscTryMethod(J,"MatSNESMFSetReuseBase_C",(Mat,PetscBool),(J,use));
156:   return(0);
157: }

159: static PetscErrorCode  MatSNESMFGetReuseBase_SNESMF(Mat J,PetscBool *use)
160: {
162:   if (J->ops->assemblyend == MatAssemblyEnd_SNESMF_UseBase) *use = PETSC_TRUE;
163:   else *use = PETSC_FALSE;
164:   return(0);
165: }

167: /*@
168:     MatSNESMFGetReuseBase - Determines if the base vector is to be used for differencing even if the function provided to SNESSetFunction() is not the
169:                        same as that provided to MatMFFDSetFunction().

171:     Logically Collective on Mat

173:     Input Parameter:
174: .   J - the MatMFFD matrix

176:     Output Parameter:
177: .   use - if true always reuse the base vector instead of recomputing f(u) even if the function in the MatSNESMF is 
178:           not SNESComputeFunction()

180:     Notes: Care must be taken when using this routine to insure that the function provided to MatMFFDSetFunction(), call it F_MF() is compatible with
181:            with that provided to SNESSetFunction(), call it F_SNES(). That is, (F_MF(u + h*d) - F_SNES(u))/h has to approximate the derivative

183:     Developer Notes: This was provided for the MOOSE team who desired to have a SNESSetFunction() function that could change configurations due
184:                      to contacts while the function provided to MatMFFDSetFunction() cannot. Except for the possibility of changing the configuration
185:                      both functions compute the same mathematical function so the differencing makes sense.

187:     Level: advanced

189: .seealso: MatCreateSNESMF(), MatSNESMFSetReuseBase()
190: @*/
191: PetscErrorCode  MatSNESMFGetReuseBase(Mat J,PetscBool *use)
192: {

197:   PetscUseMethod(J,"MatSNESMFGetReuseBase_C",(Mat,PetscBool*),(J,use));
198:   return(0);
199: }

201: /*@
202:    MatCreateSNESMF - Creates a matrix-free matrix context for use with
203:    a SNES solver.  This matrix can be used as the Jacobian argument for
204:    the routine SNESSetJacobian(). See MatCreateMFFD() for details on how
205:    the finite difference computation is done.

207:    Collective on SNES and Vec

209:    Input Parameters:
210: .  snes - the SNES context

212:    Output Parameter:
213: .  J - the matrix-free matrix

215:    Level: advanced


218:    Notes:
219:      You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
220:      preconditioner matrix

222:      If you wish to provide a different function to do differencing on to compute the matrix-free operator than
223:      that provided to SNESSetFunction() then call MatMFFDSetFunction() with your function after this call.

225:      The difference between this routine and MatCreateMFFD() is that this matrix
226:      automatically gets the current base vector from the SNES object and not from an
227:      explicit call to MatMFFDSetBase().

229:    Warning:
230:      If MatMFFDSetBase() is ever called on jac then this routine will NO longer get
231:      the x from the SNES object and MatMFFDSetBase() must from that point on be used to
232:      change the base vector x.

234:    Warning:
235:      Using a different function for the differencing will not work if you are using non-linear left preconditioning.


238: .seealso: MatDestroy(), MatMFFDSetFunction(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin()
239:           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateMFFD(),
240:           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian(), MatSNESMFSetReuseBase(), MatSNESMFGetReuseBase()

242: @*/
243: PetscErrorCode  MatCreateSNESMF(SNES snes,Mat *J)
244: {
246:   PetscInt       n,N;
247:   MatMFFD        mf;

250:   if (snes->vec_func) {
251:     VecGetLocalSize(snes->vec_func,&n);
252:     VecGetSize(snes->vec_func,&N);
253:   } else if (snes->dm) {
254:     Vec tmp;
255:     DMGetGlobalVector(snes->dm,&tmp);
256:     VecGetLocalSize(tmp,&n);
257:     VecGetSize(tmp,&N);
258:     DMRestoreGlobalVector(snes->dm,&tmp);
259:   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
260:   MatCreateMFFD(PetscObjectComm((PetscObject)snes),n,n,N,N,J);

262:   mf      = (MatMFFD)(*J)->data;
263:   mf->ctx = snes;

265:   if (snes->npc && snes->npcside== PC_LEFT) {
266:     MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunctionDefaultNPC,snes);
267:   } else {
268:     MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
269:   }

271:   (*J)->ops->assemblyend = MatAssemblyEnd_SNESMF;

273:   PetscObjectComposeFunction((PetscObject)*J,"MatMFFDSetBase_C",MatMFFDSetBase_SNESMF);
274:   PetscObjectComposeFunction((PetscObject)*J,"MatSNESMFSetReuseBase_C",MatSNESMFSetReuseBase_SNESMF);
275:   PetscObjectComposeFunction((PetscObject)*J,"MatSNESMFGetReuseBase_C",MatSNESMFGetReuseBase_SNESMF);
276:   return(0);
277: }