Actual source code: snesj2.c
petsc-3.6.1 2015-08-06
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)
22: Level: intermediate
24: Options Database Key:
25: + -snes_fd_color_use_mat - use a matrix coloring from the explicit matrix nonzero pattern instead of from the DM providing the matrix
26: . -snes_fd_color - Activates SNESComputeJacobianDefaultColor() in SNESSetFromOptions()
27: . -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
28: . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
29: - -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS)
31: Notes: If the coloring is not provided through the context, this will first try to get the
32: coloring from the DM. If the DM type has no coloring routine, then it will try to
33: get the coloring from the matrix. This requires that the matrix have nonzero entries
34: precomputed. This is discouraged, as MatColoringApply() is not parallel by default.
36: .keywords: SNES, finite differences, Jacobian, coloring, sparse
38: .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault()
39: MatFDColoringCreate(), MatFDColoringSetFunction()
41: @*/
43: PetscErrorCode SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
44: {
45: MatFDColoring color = (MatFDColoring)ctx;
47: DM dm;
48: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
49: Vec F;
50: void *funcctx;
51: MatColoring mc;
52: ISColoring iscoloring;
53: PetscBool hascolor;
54: PetscBool solvec,matcolor = PETSC_FALSE;
58: else {PetscObjectQuery((PetscObject)B,"SNESMatFDColoring",(PetscObject*)&color);}
59: SNESGetFunction(snes,&F,&func,&funcctx);
60: if (!color) {
61: SNESGetDM(snes,&dm);
62: DMHasColoring(dm,&hascolor);
63: matcolor = PETSC_FALSE;
64: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);
65: if (hascolor && !matcolor) {
66: DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);
67: MatFDColoringCreate(B,iscoloring,&color);
68: MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);
69: MatFDColoringSetFromOptions(color);
70: MatFDColoringSetUp(B,iscoloring,color);
71: ISColoringDestroy(&iscoloring);
72: } else {
73: MatColoringCreate(B,&mc);
74: MatColoringSetDistance(mc,2);
75: MatColoringSetType(mc,MATCOLORINGSL);
76: MatColoringSetFromOptions(mc);
77: MatColoringApply(mc,&iscoloring);
78: MatColoringDestroy(&mc);
79: MatFDColoringCreate(B,iscoloring,&color);
80: MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);
81: MatFDColoringSetFromOptions(color);
82: MatFDColoringSetUp(B,iscoloring,color);
83: ISColoringDestroy(&iscoloring);
84: }
85: PetscObjectCompose((PetscObject)B,"SNESMatFDColoring",(PetscObject)color);
86: PetscObjectDereference((PetscObject)color);
87: }
89: /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */
90: VecEqual(x1,snes->vec_sol,&solvec);
91: if (!snes->vec_rhs && solvec) {
92: MatFDColoringSetF(color,F);
93: }
94: MatFDColoringApply(B,color,x1,snes);
95: if (J != B) {
96: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
97: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
98: }
99: return(0);
100: }