Actual source code: snesj2.c
1: #define PETSCSNES_DLL
3: #include private/snesimpl.h
7: /*@C
8: SNESDefaultComputeJacobianColor - Computes the Jacobian using
9: finite differences and coloring to exploit matrix sparsity.
10:
11: Collective on SNES
13: Input Parameters:
14: + snes - nonlinear solver object
15: . x1 - location at which to evaluate Jacobian
16: - ctx - coloring context, where ctx must have type MatFDColoring,
17: as created via MatFDColoringCreate()
19: Output Parameters:
20: + J - Jacobian matrix (not altered in this routine)
21: . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
22: - flag - flag indicating whether the matrix sparsity structure has changed
24: Level: intermediate
26: .keywords: SNES, finite differences, Jacobian, coloring, sparse
28: .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian()
29: TSDefaultComputeJacobianColor(), MatFDColoringCreate(),
30: MatFDColoringSetFunction()
32: @*/
34: PetscErrorCode SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
35: {
36: MatFDColoring color = (MatFDColoring) ctx;
38: Vec f;
39: PetscErrorCode (*ff)(void),(*fd)(void);
42: *flag = SAME_NONZERO_PATTERN;
43: SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);
44: MatFDColoringGetFunction(color,&fd,PETSC_NULL);
45: if (fd == ff) { /* reuse function value computed in SNES */
46: MatFDColoringSetF(color,f);
47: }
48: MatFDColoringApply(*B,color,x1,flag,snes);
49: if (*J != *B) {
50: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
51: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
52: }
53: return(0);
54: }