xref: /petsc/src/snes/interface/snesj2.c (revision 3e1910f1ab6113d8365e15c6b8c907ccce7ce4ea)
1 
2 #include <petsc-private/snesimpl.h>    /*I  "petscsnes.h"  I*/
3 #include <petscdm.h>                   /*I  "petscdm.h"    I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "SNESComputeJacobianDefaultColor"
7 /*@C
8     SNESComputeJacobianDefaultColor - 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 - ignored context parameter
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 .notes: This will first try to get the coloring from the DM.  If the DM type
26         has no coloring routine, then it will try to get the coloring from the matrix.  This
27         requires that the matrix have nonzero entries precomputed, such as in
28         snes/examples/tutorials/ex45.c.
29 
30 .keywords: SNES, finite differences, Jacobian, coloring, sparse
31 
32 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault()
33           MatFDColoringCreate(), MatFDColoringSetFunction()
34 
35 @*/
36 
37 PetscErrorCode  SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
38 {
39   MatFDColoring  color = NULL;
40   PetscErrorCode ierr;
41   DM             dm;
42   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
43   Vec            F;
44   void           *funcctx;
45   ISColoring     iscoloring;
46   PetscBool      hascolor;
47 
48   PetscFunctionBegin;
49   ierr  = PetscObjectQuery((PetscObject)*B,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr);
50   *flag = SAME_NONZERO_PATTERN;
51   ierr  = SNESGetFunction(snes,&F,&func,&funcctx);CHKERRQ(ierr);
52   if (!color) {
53     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
54     ierr = DMHasColoring(dm,&hascolor);CHKERRQ(ierr);
55     if (hascolor) {
56       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);CHKERRQ(ierr);
57       ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
58       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
59       ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);CHKERRQ(ierr);
60       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
61     } else {
62       ierr = MatGetColoring(*B,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr);
63       ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
64       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
65       ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);CHKERRQ(ierr);
66       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
67     }
68     ierr = PetscObjectCompose((PetscObject)*B,"SNESMatFDColoring",(PetscObject)color);CHKERRQ(ierr);
69     ierr = PetscObjectDereference((PetscObject)color);CHKERRQ(ierr);
70   }
71 
72   /* F is only usable if there is no RHS on the SNES */
73   if (!snes->vec_rhs) {
74     ierr = MatFDColoringSetF(color,F);CHKERRQ(ierr);
75   }
76   ierr = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
77   if (*J != *B) {
78     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80   }
81   PetscFunctionReturn(0);
82 }
83