xref: /petsc/src/snes/interface/snesj.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
8*f6dfbefdSBarry Smith    Collective on snes
9fee21e36SBarry Smith 
10c7afd0dbSLois Curfman McInnes    Input Parameters:
11*f6dfbefdSBarry Smith +  snes - the `SNES` context
126b867d5aSJose E. Roman .  x1 - compute Jacobian at this point
13*f6dfbefdSBarry Smith -  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 
19*f6dfbefdSBarry Smith    Options Database Keys:
20*f6dfbefdSBarry 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)
23*f6dfbefdSBarry 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
28*f6dfbefdSBarry 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
33*f6dfbefdSBarry Smith    `SNESComputeJacobianDefaultColor()`.
34b4fc646aSLois Curfman McInnes 
35*f6dfbefdSBarry Smith    This routine ignores the maximum number of function evaluations set with `SNESSetTolerances()` and the function
36*f6dfbefdSBarry Smith    evaluations it performs are not counted in what is returned by of `SNESGetNumberFunctionEvals()`.
37*f6dfbefdSBarry Smith 
38*f6dfbefdSBarry Smith    This function can be provided to `SNESSetJacobian()` along with a dense matrix to hold the Jacobian
390df40c35SBarry Smith 
4036851e7fSLois Curfman McInnes    Level: intermediate
4136851e7fSLois Curfman McInnes 
42*f6dfbefdSBarry Smith .seealso: `SNES`, `SNESSetJacobian()`, `SNESSetJacobian()`, `SNESComputeJacobianDefaultColor()`, `MatCreateSNESMF()`
4311320018SBarry Smith @*/
449371c9d4SSatish Balay PetscErrorCode SNESComputeJacobianDefault(SNES snes, Vec x1, Mat J, Mat B, void *ctx) {
4588c956adSLois Curfman McInnes   Vec                j1a, j2a, x2;
460df40c35SBarry Smith   PetscInt           i, N, start, end, j, value, root, max_funcs = snes->max_funcs;
475edff71fSBarry Smith   PetscScalar        dx, *y, wscale;
485edff71fSBarry Smith   const PetscScalar *xx;
4977d8c4bbSBarry Smith   PetscReal          amax, epsilon = PETSC_SQRT_MACHINE_EPSILON;
5079f36460SBarry Smith   PetscReal          dx_min = 1.e-16, dx_par = 1.e-1, unorm;
51bbb6d6a8SBarry Smith   MPI_Comm           comm;
52ace3abfcSBarry Smith   PetscBool          assembled, use_wp = PETSC_TRUE, flg;
5379f36460SBarry Smith   const char        *list[2] = {"ds", "wp"};
5486c88abdSHong Zhang   PetscMPIInt        size;
5586c88abdSHong Zhang   const PetscInt    *ranges;
560df40c35SBarry Smith   DM                 dm;
570df40c35SBarry Smith   DMSNES             dms;
580521c3abSLois Curfman McInnes 
593a40ed3dSBarry Smith   PetscFunctionBegin;
600df40c35SBarry Smith   snes->max_funcs = PETSC_MAX_INT;
615d5a84c9SBarry Smith   /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */
629566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_test_err", &epsilon, NULL));
6423242f5aSBarry Smith 
659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)x1, &comm));
669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
679566063dSJacob Faibussowitsch   PetscCall(MatAssembled(B, &assembled));
681baa6e33SBarry Smith   if (assembled) PetscCall(MatZeroEntries(B));
69aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
70b1f624c7SBarry Smith     if (snes->dm) {
719566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(snes->dm, &j1a));
729566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(snes->dm, &j2a));
739566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(snes->dm, &x2));
74b1f624c7SBarry Smith     } else {
75aa79bc6dSLois Curfman McInnes       snes->nvwork = 3;
769566063dSJacob Faibussowitsch       PetscCall(VecDuplicateVecs(x1, snes->nvwork, &snes->vwork));
779566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParents(snes, snes->nvwork, snes->vwork));
789371c9d4SSatish Balay       j1a = snes->vwork[0];
799371c9d4SSatish Balay       j2a = snes->vwork[1];
809371c9d4SSatish Balay       x2  = snes->vwork[2];
81b1f624c7SBarry Smith     }
82b1f624c7SBarry Smith   }
8323242f5aSBarry Smith 
849566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x1, &N));
859566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x1, &start, &end));
869566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
879566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &dms));
880df40c35SBarry Smith   if (dms->ops->computemffunction) {
899566063dSJacob Faibussowitsch     PetscCall(SNESComputeMFFunction(snes, x1, j1a));
900df40c35SBarry Smith   } else {
919566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, x1, j1a));
920df40c35SBarry Smith   }
93c005e166SLois Curfman McInnes 
94d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing options", "SNES");
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_fd_type", "Algorithm to compute difference parameter", "SNESComputeJacobianDefault", list, 2, "wp", &value, &flg));
96d0609cedSBarry Smith   PetscOptionsEnd();
97f5af7f23SKarl Rupp   if (flg && !value) use_wp = PETSC_FALSE;
98f5af7f23SKarl Rupp 
9948a46eb9SPierre Jolivet   if (use_wp) PetscCall(VecNorm(x1, NORM_2, &unorm));
100c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
10188c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
10288c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
103c005e166SLois Curfman McInnes    */
10439e2f89bSBarry Smith   for (i = 0; i < N; i++) {
1059566063dSJacob Faibussowitsch     PetscCall(VecCopy(x1, x2));
10623242f5aSBarry Smith     if (i >= start && i < end) {
1079566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(x1, &xx));
1085620d6dcSBarry Smith       if (use_wp) dx = PetscSqrtReal(1.0 + unorm);
109f5af7f23SKarl Rupp       else dx = xx[i - start];
1109566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(x1, &xx));
1116bd79f97SJed Brown       if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par;
11239e2f89bSBarry Smith       dx *= epsilon;
11374f6f00dSLois Curfman McInnes       wscale = 1.0 / dx;
1149566063dSJacob Faibussowitsch       PetscCall(VecSetValues(x2, 1, &i, &dx, ADD_VALUES));
1156c783aadSBarry Smith     } else {
116bbb6d6a8SBarry Smith       wscale = 0.0;
117bbb6d6a8SBarry Smith     }
1189566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(x2));
1199566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(x2));
1200df40c35SBarry Smith     if (dms->ops->computemffunction) {
1219566063dSJacob Faibussowitsch       PetscCall(SNESComputeMFFunction(snes, x2, j2a));
1220df40c35SBarry Smith     } else {
1239566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, x2, j2a));
1240df40c35SBarry Smith     }
1259566063dSJacob Faibussowitsch     PetscCall(VecAXPY(j2a, -1.0, j1a));
12686c88abdSHong Zhang     /* Communicate scale=1/dx_i to all processors */
1279566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(x1, &ranges));
12886c88abdSHong Zhang     root = size;
12986c88abdSHong Zhang     for (j = size - 1; j > -1; j--) {
13086c88abdSHong Zhang       root--;
13186c88abdSHong Zhang       if (i >= ranges[j]) break;
13286c88abdSHong Zhang     }
1339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&wscale, 1, MPIU_SCALAR, root, comm));
1349566063dSJacob Faibussowitsch     PetscCall(VecScale(j2a, wscale));
1359371c9d4SSatish Balay     PetscCall(VecNorm(j2a, NORM_INFINITY, &amax));
1369371c9d4SSatish Balay     amax *= 1.e-14;
1379566063dSJacob Faibussowitsch     PetscCall(VecGetArray(j2a, &y));
13823242f5aSBarry Smith     for (j = start; j < end; j++) {
13948a46eb9SPierre Jolivet       if (PetscAbsScalar(y[j - start]) > amax || j == i) PetscCall(MatSetValues(B, 1, &j, 1, &i, y + j - start, INSERT_VALUES));
14023242f5aSBarry Smith     }
1419566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(j2a, &y));
14223242f5aSBarry Smith   }
143b1f624c7SBarry Smith   if (snes->dm) {
1449566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(snes->dm, &j1a));
1459566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(snes->dm, &j2a));
1469566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(snes->dm, &x2));
147b1f624c7SBarry Smith   }
1489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
15094ab13aaSBarry Smith   if (B != J) {
1519566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
153f588057bSBarry Smith   }
1540df40c35SBarry Smith   snes->max_funcs = max_funcs;
1550df40c35SBarry Smith   snes->nfuncs -= N;
1563a40ed3dSBarry Smith   PetscFunctionReturn(0);
15711320018SBarry Smith }
158