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