xref: /petsc/src/snes/interface/snesj2.c (revision db31f6defe33e3dfd9985704e6b8f54a7ffb6553)
1 
2 #include <petsc-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 = PETSC_NULL;
35   PetscErrorCode ierr;
36   DM dm;
37   PetscErrorCode              (*func)(SNES,Vec,Vec,void*);
38   Vec                         F;
39   void                        *funcctx;
40   ISColoring                  iscoloring;
41 
42   PetscFunctionBegin;
43   /* PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,6);*/
44   ierr = PetscObjectQuery((PetscObject)*B,"SNESMatFDColoring",(PetscObject *)&color);CHKERRQ(ierr);
45   *flag = SAME_NONZERO_PATTERN;
46   ierr = SNESGetFunction(snes,&F,&func,&funcctx);
47   if (!color) {
48     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
49     /* temporarily reset the context by creating a new one*/
50     ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);CHKERRQ(ierr);
51     ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
52     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
53     ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);CHKERRQ(ierr);
54     ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
55     ierr = PetscObjectCompose((PetscObject)*B,"SNESMatFDColoring",(PetscObject)color);CHKERRQ(ierr);
56   }
57   ierr  = MatFDColoringSetF(color,PETSC_NULL);CHKERRQ(ierr);
58   ierr  = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
59   if (*J != *B) {
60     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62   }
63   PetscFunctionReturn(0);
64 }
65