xref: /petsc/src/snes/interface/snesj.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 
2 #include <petsc/private/snesimpl.h> /*I  "petscsnes.h"  I*/
3 #include <petscdm.h>
4 
5 /*@C
6    SNESComputeJacobianDefault - Computes the Jacobian using finite differences.
7 
8    Collective on snes
9 
10    Input Parameters:
11 +  snes - the `SNES` context
12 .  x1 - compute Jacobian at this point
13 -  ctx - application's function context, as set with `SNESSetFunction()`
14 
15    Output Parameters:
16 +  J - Jacobian matrix (not altered in this routine)
17 -  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
18 
19    Options Database Keys:
20 +  -snes_fd - Activates `SNESComputeJacobianDefault()`
21 .  -snes_test_err - Square root of function error tolerance, default square root of machine
22                     epsilon (1.e-8 in double, 3.e-4 in single)
23 -  -mat_fd_type - Either wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
24 
25    Notes:
26    This routine is slow and expensive, and is not currently optimized
27    to take advantage of sparsity in the problem.  Although
28    `SNESComputeJacobianDefault()` is not recommended for general use
29    in large-scale applications, It can be useful in checking the
30    correctness of a user-provided Jacobian.
31 
32    An alternative routine that uses coloring to exploit matrix sparsity is
33    `SNESComputeJacobianDefaultColor()`.
34 
35    This routine ignores the maximum number of function evaluations set with `SNESSetTolerances()` and the function
36    evaluations it performs are not counted in what is returned by of `SNESGetNumberFunctionEvals()`.
37 
38    This function can be provided to `SNESSetJacobian()` along with a dense matrix to hold the Jacobian
39 
40    Level: intermediate
41 
42 .seealso: `SNES`, `SNESSetJacobian()`, `SNESSetJacobian()`, `SNESComputeJacobianDefaultColor()`, `MatCreateSNESMF()`
43 @*/
44 PetscErrorCode SNESComputeJacobianDefault(SNES snes, Vec x1, Mat J, Mat B, void *ctx) {
45   Vec                j1a, j2a, x2;
46   PetscInt           i, N, start, end, j, value, root, max_funcs = snes->max_funcs;
47   PetscScalar        dx, *y, wscale;
48   const PetscScalar *xx;
49   PetscReal          amax, epsilon = PETSC_SQRT_MACHINE_EPSILON;
50   PetscReal          dx_min = 1.e-16, dx_par = 1.e-1, unorm;
51   MPI_Comm           comm;
52   PetscBool          assembled, use_wp = PETSC_TRUE, flg;
53   const char        *list[2] = {"ds", "wp"};
54   PetscMPIInt        size;
55   const PetscInt    *ranges;
56   DM                 dm;
57   DMSNES             dms;
58 
59   PetscFunctionBegin;
60   snes->max_funcs = PETSC_MAX_INT;
61   /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */
62   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
63   PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_test_err", &epsilon, NULL));
64 
65   PetscCall(PetscObjectGetComm((PetscObject)x1, &comm));
66   PetscCallMPI(MPI_Comm_size(comm, &size));
67   PetscCall(MatAssembled(B, &assembled));
68   if (assembled) PetscCall(MatZeroEntries(B));
69   if (!snes->nvwork) {
70     if (snes->dm) {
71       PetscCall(DMGetGlobalVector(snes->dm, &j1a));
72       PetscCall(DMGetGlobalVector(snes->dm, &j2a));
73       PetscCall(DMGetGlobalVector(snes->dm, &x2));
74     } else {
75       snes->nvwork = 3;
76       PetscCall(VecDuplicateVecs(x1, snes->nvwork, &snes->vwork));
77       j1a = snes->vwork[0];
78       j2a = snes->vwork[1];
79       x2  = snes->vwork[2];
80     }
81   }
82 
83   PetscCall(VecGetSize(x1, &N));
84   PetscCall(VecGetOwnershipRange(x1, &start, &end));
85   PetscCall(SNESGetDM(snes, &dm));
86   PetscCall(DMGetDMSNES(dm, &dms));
87   if (dms->ops->computemffunction) {
88     PetscCall(SNESComputeMFFunction(snes, x1, j1a));
89   } else {
90     PetscCall(SNESComputeFunction(snes, x1, j1a));
91   }
92 
93   PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing options", "SNES");
94   PetscCall(PetscOptionsEList("-mat_fd_type", "Algorithm to compute difference parameter", "SNESComputeJacobianDefault", list, 2, "wp", &value, &flg));
95   PetscOptionsEnd();
96   if (flg && !value) use_wp = PETSC_FALSE;
97 
98   if (use_wp) PetscCall(VecNorm(x1, NORM_2, &unorm));
99   /* Compute Jacobian approximation, 1 column at a time.
100       x1 = current iterate, j1a = F(x1)
101       x2 = perturbed iterate, j2a = F(x2)
102    */
103   for (i = 0; i < N; i++) {
104     PetscCall(VecCopy(x1, x2));
105     if (i >= start && i < end) {
106       PetscCall(VecGetArrayRead(x1, &xx));
107       if (use_wp) dx = PetscSqrtReal(1.0 + unorm);
108       else dx = xx[i - start];
109       PetscCall(VecRestoreArrayRead(x1, &xx));
110       if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par;
111       dx *= epsilon;
112       wscale = 1.0 / dx;
113       PetscCall(VecSetValues(x2, 1, &i, &dx, ADD_VALUES));
114     } else {
115       wscale = 0.0;
116     }
117     PetscCall(VecAssemblyBegin(x2));
118     PetscCall(VecAssemblyEnd(x2));
119     if (dms->ops->computemffunction) {
120       PetscCall(SNESComputeMFFunction(snes, x2, j2a));
121     } else {
122       PetscCall(SNESComputeFunction(snes, x2, j2a));
123     }
124     PetscCall(VecAXPY(j2a, -1.0, j1a));
125     /* Communicate scale=1/dx_i to all processors */
126     PetscCall(VecGetOwnershipRanges(x1, &ranges));
127     root = size;
128     for (j = size - 1; j > -1; j--) {
129       root--;
130       if (i >= ranges[j]) break;
131     }
132     PetscCallMPI(MPI_Bcast(&wscale, 1, MPIU_SCALAR, root, comm));
133     PetscCall(VecScale(j2a, wscale));
134     PetscCall(VecNorm(j2a, NORM_INFINITY, &amax));
135     amax *= 1.e-14;
136     PetscCall(VecGetArray(j2a, &y));
137     for (j = start; j < end; j++) {
138       if (PetscAbsScalar(y[j - start]) > amax || j == i) PetscCall(MatSetValues(B, 1, &j, 1, &i, y + j - start, INSERT_VALUES));
139     }
140     PetscCall(VecRestoreArray(j2a, &y));
141   }
142   if (snes->dm) {
143     PetscCall(DMRestoreGlobalVector(snes->dm, &j1a));
144     PetscCall(DMRestoreGlobalVector(snes->dm, &j2a));
145     PetscCall(DMRestoreGlobalVector(snes->dm, &x2));
146   }
147   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
148   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
149   if (B != J) {
150     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
151     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
152   }
153   snes->max_funcs = max_funcs;
154   snes->nfuncs -= N;
155   PetscFunctionReturn(0);
156 }
157