xref: /petsc/src/snes/interface/snesj2.c (revision 9e7d6b0a155aa3071c171c8ef5a75198ca92b5e7)
1 #define PETSCSNES_DLL
2 
3 #include "include/private/matimpl.h"      /*I  "petscmat.h"  I*/
4 #include "include/private/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_lag_jacobian <freq> - -2 compute Jacobian at next call, then do not compute again, -1 do not compute Jacobian ever,
27          1 compute Jacobian every Newton iteration, 2 compute Jacobian every second Newton iteration etc.
28 
29     Level: intermediate
30 
31 .keywords: SNES, finite differences, Jacobian, coloring, sparse
32 
33 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian()
34           TSDefaultComputeJacobianColor(), MatFDColoringCreate(),
35           MatFDColoringSetFunction()
36 
37 @*/
38 
39 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
40 {
41   MatFDColoring  color = (MatFDColoring) ctx;
42   PetscErrorCode ierr;
43   PetscInt       freq,it;
44   Vec            f;
45   PetscErrorCode (*ff)(void),(*fd)(void);
46 
47   PetscFunctionBegin;
48   ierr = MatFDColoringGetLagJacobian(color,&freq);CHKERRQ(ierr);
49   ierr = SNESGetIterationNumber(snes,&it);CHKERRQ(ierr);
50 
51   if (freq == -1) {
52     ierr = PetscInfo(color,"Skipping Jacobian recomputation because freq is -1\n");CHKERRQ(ierr);
53     *flag = SAME_PRECONDITIONER;
54   } else if ((freq > 1) && ((it % freq))) {
55     ierr = PetscInfo2(color,"Skipping Jacobian recomputation, it %D, freq %D\n",it,freq);CHKERRQ(ierr);
56     *flag = SAME_PRECONDITIONER;
57   } else {
58     ierr = PetscInfo2(color,"Computing Jacobian, it %D, freq %D\n",it,freq);CHKERRQ(ierr);
59     *flag = SAME_NONZERO_PATTERN;
60     ierr  = SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);CHKERRQ(ierr);
61     ierr  = MatFDColoringGetFunction(color,&fd,PETSC_NULL);CHKERRQ(ierr);
62     if (fd == ff) { /* reuse function value computed in SNES */
63       ierr  = MatFDColoringSetF(color,f);CHKERRQ(ierr);
64     }
65     ierr  = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
66     if (freq == -2) {
67       ierr = PetscInfo(color,"Compute Jacobian freq was -2, converting to -1\n");CHKERRQ(ierr);
68       ierr = MatFDColoringSetLagJacobian(color,-1);CHKERRQ(ierr);
69     }
70   }
71   if (*J != *B) {
72     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74   }
75   PetscFunctionReturn(0);
76 }
77 
78 
79 
80