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