Actual source code: snesj2.c

petsc-3.5.1 2014-08-06
Report Typos and Errors
  2: #include <petsc-private/snesimpl.h>    /*I  "petscsnes.h"  I*/
  3: #include <petscdm.h>                   /*I  "petscdm.h"    I*/

  7: /*@C
  8:     SNESComputeJacobianDefaultColor - Computes the Jacobian using
  9:     finite differences and coloring to exploit matrix sparsity.

 11:     Collective on SNES

 13:     Input Parameters:
 14: +   snes - nonlinear solver object
 15: .   x1 - location at which to evaluate Jacobian
 16: -   ctx - MatFDColoring context or NULL

 18:     Output Parameters:
 19: +   J - Jacobian matrix (not altered in this routine)
 20: .   B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
 21: -   flag - flag indicating whether the matrix sparsity structure has changed

 23:     Level: intermediate

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

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

 32: .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault()
 33:           MatFDColoringCreate(), MatFDColoringSetFunction()

 35: @*/

 37: PetscErrorCode  SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
 38: {
 39:   MatFDColoring  color = (MatFDColoring)ctx;
 41:   DM             dm;
 42:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
 43:   Vec            F;
 44:   void           *funcctx;
 45:   MatColoring    mc;
 46:   ISColoring     iscoloring;
 47:   PetscBool      hascolor;
 48:   PetscBool      solvec,matcolor = PETSC_FALSE;

 52:   else {PetscObjectQuery((PetscObject)B,"SNESMatFDColoring",(PetscObject*)&color);}
 53:   SNESGetFunction(snes,&F,&func,&funcctx);
 54:   if (!color) {
 55:     SNESGetDM(snes,&dm);
 56:     DMHasColoring(dm,&hascolor);
 57:     matcolor = PETSC_FALSE;
 58:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);
 59:     if (hascolor && !matcolor) {
 60:       DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);
 61:       MatFDColoringCreate(B,iscoloring,&color);
 62:       MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);
 63:       MatFDColoringSetFromOptions(color);
 64:       MatFDColoringSetUp(B,iscoloring,color);
 65:       ISColoringDestroy(&iscoloring);
 66:     } else {
 67:       MatColoringCreate(B,&mc);
 68:       MatColoringSetDistance(mc,2);
 69:       MatColoringSetType(mc,MATCOLORINGSL);
 70:       MatColoringSetFromOptions(mc);
 71:       MatColoringApply(mc,&iscoloring);
 72:       MatColoringDestroy(&mc);
 73:       MatFDColoringCreate(B,iscoloring,&color);
 74:       MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);
 75:       MatFDColoringSetFromOptions(color);
 76:       MatFDColoringSetUp(B,iscoloring,color);
 77:       ISColoringDestroy(&iscoloring);
 78:     }
 79:     PetscObjectCompose((PetscObject)B,"SNESMatFDColoring",(PetscObject)color);
 80:     PetscObjectDereference((PetscObject)color);
 81:   }

 83:   /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */
 84:   VecEqual(x1,snes->vec_sol,&solvec);
 85:   if (!snes->vec_rhs && solvec) {
 86:     MatFDColoringSetF(color,F);
 87:   }
 88:   MatFDColoringApply(B,color,x1,snes);
 89:   if (J != B) {
 90:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
 91:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
 92:   }
 93:   return(0);
 94: }