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 { 46 Vec j1a, j2a, x2; 47 PetscInt i, N, start, end, j, value, root, max_funcs = snes->max_funcs; 48 PetscScalar dx, *y, wscale; 49 const PetscScalar *xx; 50 PetscReal amax, epsilon = PETSC_SQRT_MACHINE_EPSILON; 51 PetscReal dx_min = 1.e-16, dx_par = 1.e-1, unorm; 52 MPI_Comm comm; 53 PetscBool assembled, use_wp = PETSC_TRUE, flg; 54 const char *list[2] = {"ds", "wp"}; 55 PetscMPIInt size; 56 const PetscInt *ranges; 57 DM dm; 58 DMSNES dms; 59 60 PetscFunctionBegin; 61 snes->max_funcs = PETSC_MAX_INT; 62 /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 63 PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 64 PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_test_err", &epsilon, NULL)); 65 66 PetscCall(PetscObjectGetComm((PetscObject)x1, &comm)); 67 PetscCallMPI(MPI_Comm_size(comm, &size)); 68 PetscCall(MatAssembled(B, &assembled)); 69 if (assembled) PetscCall(MatZeroEntries(B)); 70 if (!snes->nvwork) { 71 if (snes->dm) { 72 PetscCall(DMGetGlobalVector(snes->dm, &j1a)); 73 PetscCall(DMGetGlobalVector(snes->dm, &j2a)); 74 PetscCall(DMGetGlobalVector(snes->dm, &x2)); 75 } else { 76 snes->nvwork = 3; 77 PetscCall(VecDuplicateVecs(x1, snes->nvwork, &snes->vwork)); 78 j1a = snes->vwork[0]; 79 j2a = snes->vwork[1]; 80 x2 = snes->vwork[2]; 81 } 82 } 83 84 PetscCall(VecGetSize(x1, &N)); 85 PetscCall(VecGetOwnershipRange(x1, &start, &end)); 86 PetscCall(SNESGetDM(snes, &dm)); 87 PetscCall(DMGetDMSNES(dm, &dms)); 88 if (dms->ops->computemffunction) { 89 PetscCall(SNESComputeMFFunction(snes, x1, j1a)); 90 } else { 91 PetscCall(SNESComputeFunction(snes, x1, j1a)); 92 } 93 94 PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing options", "SNES"); 95 PetscCall(PetscOptionsEList("-mat_fd_type", "Algorithm to compute difference parameter", "SNESComputeJacobianDefault", list, 2, "wp", &value, &flg)); 96 PetscOptionsEnd(); 97 if (flg && !value) use_wp = PETSC_FALSE; 98 99 if (use_wp) PetscCall(VecNorm(x1, NORM_2, &unorm)); 100 /* Compute Jacobian approximation, 1 column at a time. 101 x1 = current iterate, j1a = F(x1) 102 x2 = perturbed iterate, j2a = F(x2) 103 */ 104 for (i = 0; i < N; i++) { 105 PetscCall(VecCopy(x1, x2)); 106 if (i >= start && i < end) { 107 PetscCall(VecGetArrayRead(x1, &xx)); 108 if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 109 else dx = xx[i - start]; 110 PetscCall(VecRestoreArrayRead(x1, &xx)); 111 if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 112 dx *= epsilon; 113 wscale = 1.0 / dx; 114 PetscCall(VecSetValues(x2, 1, &i, &dx, ADD_VALUES)); 115 } else { 116 wscale = 0.0; 117 } 118 PetscCall(VecAssemblyBegin(x2)); 119 PetscCall(VecAssemblyEnd(x2)); 120 if (dms->ops->computemffunction) { 121 PetscCall(SNESComputeMFFunction(snes, x2, j2a)); 122 } else { 123 PetscCall(SNESComputeFunction(snes, x2, j2a)); 124 } 125 PetscCall(VecAXPY(j2a, -1.0, j1a)); 126 /* Communicate scale=1/dx_i to all processors */ 127 PetscCall(VecGetOwnershipRanges(x1, &ranges)); 128 root = size; 129 for (j = size - 1; j > -1; j--) { 130 root--; 131 if (i >= ranges[j]) break; 132 } 133 PetscCallMPI(MPI_Bcast(&wscale, 1, MPIU_SCALAR, root, comm)); 134 PetscCall(VecScale(j2a, wscale)); 135 PetscCall(VecNorm(j2a, NORM_INFINITY, &amax)); 136 amax *= 1.e-14; 137 PetscCall(VecGetArray(j2a, &y)); 138 for (j = start; j < end; j++) { 139 if (PetscAbsScalar(y[j - start]) > amax || j == i) PetscCall(MatSetValues(B, 1, &j, 1, &i, y + j - start, INSERT_VALUES)); 140 } 141 PetscCall(VecRestoreArray(j2a, &y)); 142 } 143 if (snes->dm) { 144 PetscCall(DMRestoreGlobalVector(snes->dm, &j1a)); 145 PetscCall(DMRestoreGlobalVector(snes->dm, &j2a)); 146 PetscCall(DMRestoreGlobalVector(snes->dm, &x2)); 147 } 148 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 149 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 150 if (B != J) { 151 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 152 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 153 } 154 snes->max_funcs = max_funcs; 155 snes->nfuncs -= N; 156 PetscFunctionReturn(0); 157 } 158