Actual source code: snesj.c
1: /*$Id: snesj.c,v 1.75 2001/09/11 18:06:40 bsmith Exp $*/
3: #include src/snes/snesimpl.h
5: #undef __FUNCT__
7: /*@C
8: SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
10: Collective on SNES
12: Input Parameters:
13: + x1 - compute Jacobian at this point
14: - ctx - application's function context, as set with SNESSetFunction()
16: Output Parameters:
17: + J - Jacobian matrix (not altered in this routine)
18: . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
19: - flag - flag indicating whether the matrix sparsity structure has changed
21: Options Database Key:
22: + -snes_fd - Activates SNESDefaultComputeJacobian()
23: - -snes_test_err - Square root of function error tolerance, default square root of machine
24: epsilon (1.e-8 in double, 3.e-4 in single)
26: Notes:
27: This routine is slow and expensive, and is not currently optimized
28: to take advantage of sparsity in the problem. Although
29: SNESDefaultComputeJacobian() is not recommended for general use
30: in large-scale applications, It can be useful in checking the
31: correctness of a user-provided Jacobian.
33: An alternative routine that uses coloring to explot matrix sparsity is
34: SNESDefaultComputeJacobianColor().
36: Level: intermediate
38: .keywords: SNES, finite differences, Jacobian
40: .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor()
41: @*/
42: int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
43: {
44: Vec j1a,j2a,x2;
45: int i,ierr,N,start,end,j;
46: PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
47: PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
48: PetscReal dx_min = 1.e-16,dx_par = 1.e-1;
49: MPI_Comm comm;
50: int (*eval_fct)(SNES,Vec,Vec)=0;
53: PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);
54: eval_fct = SNESComputeFunction;
56: PetscObjectGetComm((PetscObject)x1,&comm);
57: MatZeroEntries(*B);
58: if (!snes->nvwork) {
59: VecDuplicateVecs(x1,3,&snes->vwork);
60: snes->nvwork = 3;
61: PetscLogObjectParents(snes,3,snes->vwork);
62: }
63: j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
65: VecGetSize(x1,&N);
66: VecGetOwnershipRange(x1,&start,&end);
67: (*eval_fct)(snes,x1,j1a);
69: /* Compute Jacobian approximation, 1 column at a time.
70: x1 = current iterate, j1a = F(x1)
71: x2 = perturbed iterate, j2a = F(x2)
72: */
73: for (i=0; i<N; i++) {
74: VecCopy(x1,x2);
75: if (i>= start && i<end) {
76: VecGetArray(x1,&xx);
77: dx = xx[i-start];
78: VecRestoreArray(x1,&xx);
79: #if !defined(PETSC_USE_COMPLEX)
80: if (dx < dx_min && dx >= 0.0) dx = dx_par;
81: else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
82: #else
83: if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
84: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
85: #endif
86: dx *= epsilon;
87: wscale = 1.0/dx;
88: VecSetValues(x2,1,&i,&dx,ADD_VALUES);
89: } else {
90: wscale = 0.0;
91: }
92: (*eval_fct)(snes,x2,j2a);
93: VecAXPY(&mone,j1a,j2a);
94: /* Communicate scale to all processors */
95: MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);
96: VecScale(&scale,j2a);
97: VecNorm(j2a,NORM_INFINITY,&amax); amax *= 1.e-14;
98: VecGetArray(j2a,&y);
99: for (j=start; j<end; j++) {
100: if (PetscAbsScalar(y[j-start]) > amax) {
101: MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);
102: }
103: }
104: VecRestoreArray(j2a,&y);
105: }
106: ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
107: ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
108: *flag = DIFFERENT_NONZERO_PATTERN;
109: return(0);
110: }