xref: /petsc/src/snes/interface/snesj2.c (revision f68b968ce39302dfa79eb1a6cfa1998ce074e829)
1 #define PETSCSNES_DLL
2 
3 #include "src/mat/matimpl.h"      /*I  "petscmat.h"  I*/
4 #include "src/snes/snesimpl.h"    /*I  "petscsnes.h"  I*/
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "SNESDefaultComputeJacobianColor"
8 /*@C
9     SNESDefaultComputeJacobianColor - Computes the Jacobian using
10     finite differences and coloring to exploit matrix sparsity.
11 
12     Collective on SNES
13 
14     Input Parameters:
15 +   snes - nonlinear solver object
16 .   x1 - location at which to evaluate Jacobian
17 -   ctx - coloring context, where ctx must have type MatFDColoring,
18           as created via MatFDColoringCreate()
19 
20     Output Parameters:
21 +   J - Jacobian matrix (not altered in this routine)
22 .   B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
23 -   flag - flag indicating whether the matrix sparsity structure has changed
24 
25     Options Database Keys:
26 .  -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor()
27 
28     Level: intermediate
29 
30 .keywords: SNES, finite differences, Jacobian, coloring, sparse
31 
32 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian()
33           TSDefaultComputeJacobianColor(), MatFDColoringCreate(),
34           MatFDColoringSetFunction()
35 
36 @*/
37 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
38 {
39   MatFDColoring  color = (MatFDColoring) ctx;
40   PetscErrorCode ierr;
41   PetscInt       freq,it;
42   Vec            f;
43   PetscErrorCode (*ff)(void),(*fd)(void);
44 
45   PetscFunctionBegin;
46   ierr = MatFDColoringGetFrequency(color,&freq);CHKERRQ(ierr);
47   ierr = SNESGetIterationNumber(snes,&it);CHKERRQ(ierr);
48 
49   if ((freq > 1) && ((it % freq))) {
50     ierr = PetscLogInfo((color,"SNESDefaultComputeJacobianColor:Skipping Jacobian recomputation, it %D, freq %D\n",it,freq));CHKERRQ(ierr);
51     *flag = SAME_PRECONDITIONER;
52   } else {
53     ierr = PetscLogInfo((color,"SNESDefaultComputeJacobianColor:Computing Jacobian, it %D, freq %D\n",it,freq));CHKERRQ(ierr);
54     *flag = SAME_NONZERO_PATTERN;
55     ierr  = SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);CHKERRQ(ierr);
56     ierr  = MatFDColoringGetFunction(color,&fd,PETSC_NULL);CHKERRQ(ierr);
57     if (fd == ff) { /* reuse function value computed in SNES */
58       ierr  = MatFDColoringSetF(color,f);CHKERRQ(ierr);
59     }
60     ierr  = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
61   }
62   if (*J != *B) {
63     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65   }
66   PetscFunctionReturn(0);
67 }
68 
69 
70 
71