Actual source code: snesj2.c
petsc-3.5.4 2015-05-23
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: .notes: If the coloring is not provided through the context, this will first try to get the
25: coloring from the DM. If the DM type has no coloring routine, then it will try to
26: get the coloring from the matrix. This requires that the matrix have nonzero entries
27: precomputed. This is discouraged, as MatColoringApply() is not parallel by default.
29: .keywords: SNES, finite differences, Jacobian, coloring, sparse
31: .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault()
32: MatFDColoringCreate(), MatFDColoringSetFunction()
34: @*/
36: PetscErrorCode SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
37: {
38: MatFDColoring color = (MatFDColoring)ctx;
40: DM dm;
41: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
42: Vec F;
43: void *funcctx;
44: MatColoring mc;
45: ISColoring iscoloring;
46: PetscBool hascolor;
47: PetscBool solvec,matcolor = PETSC_FALSE;
51: else {PetscObjectQuery((PetscObject)B,"SNESMatFDColoring",(PetscObject*)&color);}
52: SNESGetFunction(snes,&F,&func,&funcctx);
53: if (!color) {
54: SNESGetDM(snes,&dm);
55: DMHasColoring(dm,&hascolor);
56: matcolor = PETSC_FALSE;
57: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);
58: if (hascolor && !matcolor) {
59: DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);
60: MatFDColoringCreate(B,iscoloring,&color);
61: MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);
62: MatFDColoringSetFromOptions(color);
63: MatFDColoringSetUp(B,iscoloring,color);
64: ISColoringDestroy(&iscoloring);
65: } else {
66: MatColoringCreate(B,&mc);
67: MatColoringSetDistance(mc,2);
68: MatColoringSetType(mc,MATCOLORINGSL);
69: MatColoringSetFromOptions(mc);
70: MatColoringApply(mc,&iscoloring);
71: MatColoringDestroy(&mc);
72: MatFDColoringCreate(B,iscoloring,&color);
73: MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);
74: MatFDColoringSetFromOptions(color);
75: MatFDColoringSetUp(B,iscoloring,color);
76: ISColoringDestroy(&iscoloring);
77: }
78: PetscObjectCompose((PetscObject)B,"SNESMatFDColoring",(PetscObject)color);
79: PetscObjectDereference((PetscObject)color);
80: }
82: /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */
83: VecEqual(x1,snes->vec_sol,&solvec);
84: if (!snes->vec_rhs && solvec) {
85: MatFDColoringSetF(color,F);
86: }
87: MatFDColoringApply(B,color,x1,snes);
88: if (J != B) {
89: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
91: }
92: return(0);
93: }