xref: /petsc/src/snes/interface/snesj2.c (revision 9c30b7d2697335155d7490a7e085415ee7b4a02a)
1 
2 #include "src/mat/matimpl.h"      /*I  "petscmat.h"  I*/
3 #include "src/snes/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     Options Database Keys:
25 .  -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor()
26 
27     Level: intermediate
28 
29 .keywords: SNES, finite differences, Jacobian, coloring, sparse
30 
31 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian()
32           TSDefaultComputeJacobianColor(), MatFDColoringCreate(),
33           MatFDColoringSetFunction()
34 
35 @*/
36 int SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
37 {
38   MatFDColoring color = (MatFDColoring) ctx;
39   int           ierr,freq,it;
40   Vec           f;
41 
42   PetscFunctionBegin;
43   ierr = MatFDColoringGetFrequency(color,&freq);CHKERRQ(ierr);
44   ierr = SNESGetIterationNumber(snes,&it);CHKERRQ(ierr);
45 
46   if ((freq > 1) && ((it % freq))) {
47     PetscLogInfo(color,"SNESDefaultComputeJacobianColor:Skipping Jacobian recomputation, it %d, freq %d\n",it,freq);
48     *flag = SAME_PRECONDITIONER;
49   } else {
50     PetscLogInfo(color,"SNESDefaultComputeJacobianColor:Computing Jacobian, it %d, freq %d\n",it,freq);
51     *flag = SAME_NONZERO_PATTERN;
52     ierr  = SNESGetFunction(snes,&f,0,0);CHKERRQ(ierr);
53     ierr  = MatFDColoringSetF(color,f);CHKERRQ(ierr);
54     ierr  = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
55   }
56   if (J != B) {
57     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
58     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59   }
60   PetscFunctionReturn(0);
61 }
62 
63 
64 
65