Actual source code: snesj2.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <petsc/private/snesimpl.h>    /*I  "petscsnes.h"  I*/
  3: #include <petscdm.h>                   /*I  "petscdm.h"    I*/

  5: /*
  6:    MatFDColoringSetFunction() takes a function with four arguments, we want to use SNESComputeFunction()
  7:    since it logs function computation information.
  8: */
  9: static PetscErrorCode SNESComputeFunctionCtx(SNES snes,Vec x,Vec f,void *ctx)
 10: {
 11:   return SNESComputeFunction(snes,x,f);
 12: }

 16: /*@C
 17:     SNESComputeJacobianDefaultColor - Computes the Jacobian using
 18:     finite differences and coloring to exploit matrix sparsity.

 20:     Collective on SNES

 22:     Input Parameters:
 23: +   snes - nonlinear solver object
 24: .   x1 - location at which to evaluate Jacobian
 25: -   ctx - MatFDColoring context or NULL

 27:     Output Parameters:
 28: +   J - Jacobian matrix (not altered in this routine)
 29: -   B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)

 31:     Level: intermediate

 33:    Options Database Key:
 34: +  -snes_fd_color_use_mat - use a matrix coloring from the explicit matrix nonzero pattern instead of from the DM providing the matrix
 35: .  -snes_fd_color - Activates SNESComputeJacobianDefaultColor() in SNESSetFromOptions()
 36: .  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
 37: .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
 38: -  -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS)

 40:     Notes: If the coloring is not provided through the context, this will first try to get the
 41:         coloring from the DM.  If the DM type has no coloring routine, then it will try to
 42:         get the coloring from the matrix.  This requires that the matrix have nonzero entries
 43:         precomputed.  This is discouraged, as MatColoringApply() is not parallel by default.

 45: .keywords: SNES, finite differences, Jacobian, coloring, sparse

 47: .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault()
 48:           MatFDColoringCreate(), MatFDColoringSetFunction()

 50: @*/

 52: PetscErrorCode  SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
 53: {
 54:   MatFDColoring  color = (MatFDColoring)ctx;
 56:   DM             dm;
 57:   MatColoring    mc;
 58:   ISColoring     iscoloring;
 59:   PetscBool      hascolor;
 60:   PetscBool      solvec,matcolor = PETSC_FALSE;

 64:   else {PetscObjectQuery((PetscObject)B,"SNESMatFDColoring",(PetscObject*)&color);}

 66:   if (!color) {
 67:     SNESGetDM(snes,&dm);
 68:     DMHasColoring(dm,&hascolor);
 69:     matcolor = PETSC_FALSE;
 70:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);
 71:     if (hascolor && !matcolor) {
 72:       DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);
 73:       MatFDColoringCreate(B,iscoloring,&color);
 74:       MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))SNESComputeFunctionCtx,NULL);
 75:       MatFDColoringSetFromOptions(color);
 76:       MatFDColoringSetUp(B,iscoloring,color);
 77:       ISColoringDestroy(&iscoloring);
 78:     } else {
 79:       MatColoringCreate(B,&mc);
 80:       MatColoringSetDistance(mc,2);
 81:       MatColoringSetType(mc,MATCOLORINGSL);
 82:       MatColoringSetFromOptions(mc);
 83:       MatColoringApply(mc,&iscoloring);
 84:       MatColoringDestroy(&mc);
 85:       MatFDColoringCreate(B,iscoloring,&color);
 86:       MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))SNESComputeFunctionCtx,NULL);
 87:       MatFDColoringSetFromOptions(color);
 88:       MatFDColoringSetUp(B,iscoloring,color);
 89:       ISColoringDestroy(&iscoloring);
 90:     }
 91:     PetscObjectCompose((PetscObject)B,"SNESMatFDColoring",(PetscObject)color);
 92:     PetscObjectDereference((PetscObject)color);
 93:   }

 95:   /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */
 96:   VecEqual(x1,snes->vec_sol,&solvec);
 97:   if (!snes->vec_rhs && solvec) {
 98:     Vec F;
 99:     SNESGetFunction(snes,&F,NULL,NULL);
100:     MatFDColoringSetF(color,F);
101:   }
102:   MatFDColoringApply(B,color,x1,snes);
103:   if (J != B) {
104:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
105:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
106:   }
107:   return(0);
108: }