xref: /petsc/src/snes/interface/snesj2.c (revision f08646a875d6c8bb0dc8fb3d7c6d74e67952873a)
1 
2 #include <private/snesimpl.h>    /*I  "petscsnes.h"  I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESDefaultComputeJacobianColor"
6 /*@C
7     SNESDefaultComputeJacobianColor - Computes the Jacobian using
8     finite differences and coloring to exploit matrix sparsity.
9 
10     Collective on SNES
11 
12     Input Parameters:
13 +   snes - nonlinear solver object
14 .   x1 - location at which to evaluate Jacobian
15 -   ctx - coloring context, where ctx must have type MatFDColoring,
16           as created via MatFDColoringCreate()
17 
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
22 
23     Level: intermediate
24 
25 .keywords: SNES, finite differences, Jacobian, coloring, sparse
26 
27 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian()
28           MatFDColoringCreate(), MatFDColoringSetFunction()
29 
30 @*/
31 
32 PetscErrorCode  SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
33 {
34   MatFDColoring  color = (MatFDColoring) ctx;
35   PetscErrorCode ierr;
36   Vec            f;
37   PetscErrorCode (*ff)(void),(*fd)(void);
38 
39   PetscFunctionBegin;
40   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,6);
41   *flag = SAME_NONZERO_PATTERN;
42   ierr  = SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);CHKERRQ(ierr);
43   ierr  = MatFDColoringGetFunction(color,&fd,PETSC_NULL);CHKERRQ(ierr);
44   if (fd == ff && !snes->vec_rhs) { /* reuse function value computed in SNES */
45     ierr  = MatFDColoringSetF(color,f);CHKERRQ(ierr);
46   }
47   ierr  = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
48   if (*J != *B) {
49     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
51   }
52   PetscFunctionReturn(0);
53 }
54 
55 
56 
57