xref: /petsc/src/snes/interface/snesj.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
111320018SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I  "petscsnes.h"  I*/
3b1f624c7SBarry Smith #include <petscdm.h>
411320018SBarry Smith 
54b828684SBarry Smith /*@C
68d359177SBarry Smith    SNESComputeJacobianDefault - Computes the Jacobian using finite differences.
711320018SBarry Smith 
8fee21e36SBarry Smith    Collective on SNES
9fee21e36SBarry Smith 
10c7afd0dbSLois Curfman McInnes    Input Parameters:
116b867d5aSJose E. Roman +  snes - the SNES context
126b867d5aSJose E. Roman .  x1 - compute Jacobian at this point
13c7afd0dbSLois Curfman McInnes -  ctx - application's function context, as set with SNESSetFunction()
14c7afd0dbSLois Curfman McInnes 
15c7afd0dbSLois Curfman McInnes    Output Parameters:
16c7afd0dbSLois Curfman McInnes +  J - Jacobian matrix (not altered in this routine)
17d1e9a80fSBarry Smith -  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
18c7afd0dbSLois Curfman McInnes 
19ad960d00SLois Curfman McInnes    Options Database Key:
208d359177SBarry Smith +  -snes_fd - Activates SNESComputeJacobianDefault()
2179f36460SBarry Smith .  -snes_test_err - Square root of function error tolerance, default square root of machine
2277d8c4bbSBarry Smith                     epsilon (1.e-8 in double, 3.e-4 in single)
233ec795f1SBarry Smith -  -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS)
24ad960d00SLois Curfman McInnes 
255f3c43d9SLois Curfman McInnes    Notes:
265f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
275f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
288d359177SBarry Smith    SNESComputeJacobianDefault() is not recommended for general use
295f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
305f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3111320018SBarry Smith 
3279f36460SBarry Smith    An alternative routine that uses coloring to exploit matrix sparsity is
338d359177SBarry Smith    SNESComputeJacobianDefaultColor().
34b4fc646aSLois Curfman McInnes 
350df40c35SBarry Smith    This routine ignores the maximum number of function evaluations set with SNESSetTolerances() and the function
360df40c35SBarry Smith    evaluations it performs are not counted in what is returned by of SNESGetNumberFunctionEvals().
370df40c35SBarry Smith 
3836851e7fSLois Curfman McInnes    Level: intermediate
3936851e7fSLois Curfman McInnes 
40db781477SPatrick Sanan .seealso: `SNESSetJacobian()`, `SNESComputeJacobianDefaultColor()`, `MatCreateSNESMF()`
4111320018SBarry Smith @*/
429371c9d4SSatish Balay PetscErrorCode SNESComputeJacobianDefault(SNES snes, Vec x1, Mat J, Mat B, void *ctx) {
4388c956adSLois Curfman McInnes   Vec                j1a, j2a, x2;
440df40c35SBarry Smith   PetscInt           i, N, start, end, j, value, root, max_funcs = snes->max_funcs;
455edff71fSBarry Smith   PetscScalar        dx, *y, wscale;
465edff71fSBarry Smith   const PetscScalar *xx;
4777d8c4bbSBarry Smith   PetscReal          amax, epsilon = PETSC_SQRT_MACHINE_EPSILON;
4879f36460SBarry Smith   PetscReal          dx_min = 1.e-16, dx_par = 1.e-1, unorm;
49bbb6d6a8SBarry Smith   MPI_Comm           comm;
50ace3abfcSBarry Smith   PetscBool          assembled, use_wp = PETSC_TRUE, flg;
5179f36460SBarry Smith   const char        *list[2] = {"ds", "wp"};
5286c88abdSHong Zhang   PetscMPIInt        size;
5386c88abdSHong Zhang   const PetscInt    *ranges;
540df40c35SBarry Smith   DM                 dm;
550df40c35SBarry Smith   DMSNES             dms;
560521c3abSLois Curfman McInnes 
573a40ed3dSBarry Smith   PetscFunctionBegin;
580df40c35SBarry Smith   snes->max_funcs = PETSC_MAX_INT;
595d5a84c9SBarry Smith   /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */
609566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_test_err", &epsilon, NULL));
6223242f5aSBarry Smith 
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)x1, &comm));
649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
659566063dSJacob Faibussowitsch   PetscCall(MatAssembled(B, &assembled));
661baa6e33SBarry Smith   if (assembled) PetscCall(MatZeroEntries(B));
67aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
68b1f624c7SBarry Smith     if (snes->dm) {
699566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(snes->dm, &j1a));
709566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(snes->dm, &j2a));
719566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(snes->dm, &x2));
72b1f624c7SBarry Smith     } else {
73aa79bc6dSLois Curfman McInnes       snes->nvwork = 3;
749566063dSJacob Faibussowitsch       PetscCall(VecDuplicateVecs(x1, snes->nvwork, &snes->vwork));
759566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParents(snes, snes->nvwork, snes->vwork));
769371c9d4SSatish Balay       j1a = snes->vwork[0];
779371c9d4SSatish Balay       j2a = snes->vwork[1];
789371c9d4SSatish Balay       x2  = snes->vwork[2];
79b1f624c7SBarry Smith     }
80b1f624c7SBarry Smith   }
8123242f5aSBarry Smith 
829566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x1, &N));
839566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x1, &start, &end));
849566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
859566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &dms));
860df40c35SBarry Smith   if (dms->ops->computemffunction) {
879566063dSJacob Faibussowitsch     PetscCall(SNESComputeMFFunction(snes, x1, j1a));
880df40c35SBarry Smith   } else {
899566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, x1, j1a));
900df40c35SBarry Smith   }
91c005e166SLois Curfman McInnes 
92d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing options", "SNES");
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_fd_type", "Algorithm to compute difference parameter", "SNESComputeJacobianDefault", list, 2, "wp", &value, &flg));
94d0609cedSBarry Smith   PetscOptionsEnd();
95f5af7f23SKarl Rupp   if (flg && !value) use_wp = PETSC_FALSE;
96f5af7f23SKarl Rupp 
97*48a46eb9SPierre Jolivet   if (use_wp) PetscCall(VecNorm(x1, NORM_2, &unorm));
98c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
9988c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
10088c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
101c005e166SLois Curfman McInnes    */
10239e2f89bSBarry Smith   for (i = 0; i < N; i++) {
1039566063dSJacob Faibussowitsch     PetscCall(VecCopy(x1, x2));
10423242f5aSBarry Smith     if (i >= start && i < end) {
1059566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(x1, &xx));
1065620d6dcSBarry Smith       if (use_wp) dx = PetscSqrtReal(1.0 + unorm);
107f5af7f23SKarl Rupp       else dx = xx[i - start];
1089566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(x1, &xx));
1096bd79f97SJed Brown       if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par;
11039e2f89bSBarry Smith       dx *= epsilon;
11174f6f00dSLois Curfman McInnes       wscale = 1.0 / dx;
1129566063dSJacob Faibussowitsch       PetscCall(VecSetValues(x2, 1, &i, &dx, ADD_VALUES));
1136c783aadSBarry Smith     } else {
114bbb6d6a8SBarry Smith       wscale = 0.0;
115bbb6d6a8SBarry Smith     }
1169566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(x2));
1179566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(x2));
1180df40c35SBarry Smith     if (dms->ops->computemffunction) {
1199566063dSJacob Faibussowitsch       PetscCall(SNESComputeMFFunction(snes, x2, j2a));
1200df40c35SBarry Smith     } else {
1219566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, x2, j2a));
1220df40c35SBarry Smith     }
1239566063dSJacob Faibussowitsch     PetscCall(VecAXPY(j2a, -1.0, j1a));
12486c88abdSHong Zhang     /* Communicate scale=1/dx_i to all processors */
1259566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(x1, &ranges));
12686c88abdSHong Zhang     root = size;
12786c88abdSHong Zhang     for (j = size - 1; j > -1; j--) {
12886c88abdSHong Zhang       root--;
12986c88abdSHong Zhang       if (i >= ranges[j]) break;
13086c88abdSHong Zhang     }
1319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&wscale, 1, MPIU_SCALAR, root, comm));
1329566063dSJacob Faibussowitsch     PetscCall(VecScale(j2a, wscale));
1339371c9d4SSatish Balay     PetscCall(VecNorm(j2a, NORM_INFINITY, &amax));
1349371c9d4SSatish Balay     amax *= 1.e-14;
1359566063dSJacob Faibussowitsch     PetscCall(VecGetArray(j2a, &y));
13623242f5aSBarry Smith     for (j = start; j < end; j++) {
137*48a46eb9SPierre Jolivet       if (PetscAbsScalar(y[j - start]) > amax || j == i) PetscCall(MatSetValues(B, 1, &j, 1, &i, y + j - start, INSERT_VALUES));
13823242f5aSBarry Smith     }
1399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(j2a, &y));
14023242f5aSBarry Smith   }
141b1f624c7SBarry Smith   if (snes->dm) {
1429566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(snes->dm, &j1a));
1439566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(snes->dm, &j2a));
1449566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(snes->dm, &x2));
145b1f624c7SBarry Smith   }
1469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
14894ab13aaSBarry Smith   if (B != J) {
1499566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
151f588057bSBarry Smith   }
1520df40c35SBarry Smith   snes->max_funcs = max_funcs;
1530df40c35SBarry Smith   snes->nfuncs -= N;
1543a40ed3dSBarry Smith   PetscFunctionReturn(0);
15511320018SBarry Smith }
156