Actual source code: snesj2.c
petsc-3.4.5 2014-06-29
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 MatGetColoring() is not parallel.
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,MatStructure *flag,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: ISColoring iscoloring;
46: PetscBool hascolor;
47: PetscBool solvec;
51: else {PetscObjectQuery((PetscObject)*B,"SNESMatFDColoring",(PetscObject*)&color);}
52: *flag = SAME_NONZERO_PATTERN;
53: SNESGetFunction(snes,&F,&func,&funcctx);
54: if (!color) {
55: SNESGetDM(snes,&dm);
56: DMHasColoring(dm,&hascolor);
57: if (hascolor) {
58: DMCreateColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);
59: MatFDColoringCreate(*B,iscoloring,&color);
60: ISColoringDestroy(&iscoloring);
61: MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);
62: MatFDColoringSetFromOptions(color);
63: } else {
64: MatGetColoring(*B,MATCOLORINGSL,&iscoloring);
65: MatFDColoringCreate(*B,iscoloring,&color);
66: ISColoringDestroy(&iscoloring);
67: MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);
68: MatFDColoringSetFromOptions(color);
69: }
70: PetscObjectCompose((PetscObject)*B,"SNESMatFDColoring",(PetscObject)color);
71: PetscObjectDereference((PetscObject)color);
72: }
74: /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */
75: VecEqual(x1,snes->vec_sol,&solvec);
76: if (!snes->vec_rhs && solvec) {
77: MatFDColoringSetF(color,F);
78: }
79: MatFDColoringApply(*B,color,x1,flag,snes);
80: if (*J != *B) {
81: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
83: }
84: return(0);
85: }