Actual source code: snesj.c
2: #include <petsc/private/snesimpl.h>
3: #include <petscdm.h>
5: /*@C
6: SNESComputeJacobianDefault - Computes the Jacobian using finite differences.
8: Collective on SNES
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()
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)
19: Options Database Key:
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)
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.
32: An alternative routine that uses coloring to exploit matrix sparsity is
33: SNESComputeJacobianDefaultColor().
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().
38: Level: intermediate
40: .seealso: SNESSetJacobian(), SNESComputeJacobianDefaultColor(), MatCreateSNESMF()
41: @*/
42: PetscErrorCode SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
43: {
44: Vec j1a,j2a,x2;
45: PetscErrorCode ierr;
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;
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: MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
63: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL);
65: PetscObjectGetComm((PetscObject)x1,&comm);
66: MPI_Comm_size(comm,&size);
67: MatAssembled(B,&assembled);
68: if (assembled) {
69: MatZeroEntries(B);
70: }
71: if (!snes->nvwork) {
72: if (snes->dm) {
73: DMGetGlobalVector(snes->dm,&j1a);
74: DMGetGlobalVector(snes->dm,&j2a);
75: DMGetGlobalVector(snes->dm,&x2);
76: } else {
77: snes->nvwork = 3;
78: VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);
79: PetscLogObjectParents(snes,snes->nvwork,snes->vwork);
80: j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
81: }
82: }
84: VecGetSize(x1,&N);
85: VecGetOwnershipRange(x1,&start,&end);
86: SNESGetDM(snes,&dm);
87: DMGetDMSNES(dm,&dms);
88: if (dms->ops->computemffunction) {
89: SNESComputeMFFunction(snes,x1,j1a);
90: } else {
91: SNESComputeFunction(snes,x1,j1a);
92: }
94: PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES");
95: 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;
99: if (use_wp) {
100: VecNorm(x1,NORM_2,&unorm);
101: }
102: /* Compute Jacobian approximation, 1 column at a time.
103: x1 = current iterate, j1a = F(x1)
104: x2 = perturbed iterate, j2a = F(x2)
105: */
106: for (i=0; i<N; i++) {
107: VecCopy(x1,x2);
108: if (i>= start && i<end) {
109: VecGetArrayRead(x1,&xx);
110: if (use_wp) dx = PetscSqrtReal(1.0 + unorm);
111: else dx = xx[i-start];
112: VecRestoreArrayRead(x1,&xx);
113: if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par;
114: dx *= epsilon;
115: wscale = 1.0/dx;
116: VecSetValues(x2,1,&i,&dx,ADD_VALUES);
117: } else {
118: wscale = 0.0;
119: }
120: VecAssemblyBegin(x2);
121: VecAssemblyEnd(x2);
122: if (dms->ops->computemffunction) {
123: SNESComputeMFFunction(snes,x2,j2a);
124: } else {
125: SNESComputeFunction(snes,x2,j2a);
126: }
127: VecAXPY(j2a,-1.0,j1a);
128: /* Communicate scale=1/dx_i to all processors */
129: VecGetOwnershipRanges(x1,&ranges);
130: root = size;
131: for (j=size-1; j>-1; j--) {
132: root--;
133: if (i>=ranges[j]) break;
134: }
135: MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);
136: VecScale(j2a,wscale);
137: VecNorm(j2a,NORM_INFINITY,&amax); amax *= 1.e-14;
138: VecGetArray(j2a,&y);
139: for (j=start; j<end; j++) {
140: if (PetscAbsScalar(y[j-start]) > amax || j == i) {
141: MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES);
142: }
143: }
144: VecRestoreArray(j2a,&y);
145: }
146: if (snes->dm) {
147: DMRestoreGlobalVector(snes->dm,&j1a);
148: DMRestoreGlobalVector(snes->dm,&j2a);
149: DMRestoreGlobalVector(snes->dm,&x2);
150: }
151: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
152: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
153: if (B != J) {
154: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
155: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
156: }
157: snes->max_funcs = max_funcs;
158: snes->nfuncs -= N;
159: return(0);
160: }