xref: /petsc/src/snes/interface/snesj2.c (revision 6d75e2106b7ce71d78ffad317d596124b1823eef)
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 - ignored context parameter
16 
17     Output Parameters:
18 +   J - Jacobian matrix (not altered in this routine)
19 .   B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
20 -   flag - flag indicating whether the matrix sparsity structure has changed
21 
22     Level: intermediate
23 
24 .notes: This will first try to get the coloring from the DM.  If the DM type
25         has no coloring routine, then it will try to get the coloring from the matrix.  This
26         requires that the matrix have nonzero entries precomputed, such as in
27         snes/examples/tutorials/ex45.c.
28 
29 .keywords: SNES, finite differences, Jacobian, coloring, sparse
30 
31 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian()
32           MatFDColoringCreate(), MatFDColoringSetFunction()
33 
34 @*/
35 
36 PetscErrorCode  SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
37 {
38   MatFDColoring  color = PETSC_NULL;
39   PetscErrorCode ierr;
40   DM             dm;
41   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
42   Vec            F;
43   void           *funcctx;
44   ISColoring     iscoloring;
45   PetscBool      hascolor;
46 
47   PetscFunctionBegin;
48   ierr = PetscObjectQuery((PetscObject)*B,"SNESMatFDColoring",(PetscObject *)&color);CHKERRQ(ierr);
49   *flag = SAME_NONZERO_PATTERN;
50   ierr = SNESGetFunction(snes,&F,&func,&funcctx);CHKERRQ(ierr);
51   if (!color) {
52     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
53     ierr = DMHasColoring(dm,&hascolor);CHKERRQ(ierr);
54     if (hascolor) {
55       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);CHKERRQ(ierr);
56       ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
57       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
58       ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);CHKERRQ(ierr);
59       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
60     } else {
61       ierr = MatGetColoring(*B,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr);
62       ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
63       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
64       ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);CHKERRQ(ierr);
65       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
66     }
67     ierr = PetscObjectCompose((PetscObject)*B,"SNESMatFDColoring",(PetscObject)color);CHKERRQ(ierr);
68     ierr = PetscObjectDereference((PetscObject)color);CHKERRQ(ierr);
69   }
70 
71   /* F is only usable if there is no RHS on the SNES */
72   if (!snes->vec_rhs) {
73     ierr  = MatFDColoringSetF(color,F);CHKERRQ(ierr);
74   }
75   ierr  = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
76   if (*J != *B) {
77     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79   }
80   PetscFunctionReturn(0);
81 }
82