xref: /petsc/src/snes/interface/snesj.c (revision 9b0d118972bd5614e9e2d4aa76805104fdfe5302)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: snesj.c,v 1.13 1995/05/12 21:02:24 curfman Exp curfman $";
4 #endif
5 
6 #include "draw.h"
7 #include "snes.h"
8 
9 /*@
10    SNESDefaultComputeJacobian - Computes the Jacobian using finite
11    differences.
12 
13    Input Parameters:
14 .  x1 - compute Jacobian at this point
15 .  ctx - application's function context, as set with SNESSetFunction()
16 
17    Output Parameters:
18 .  J - Jacobian
19 .  B - preconditioner, same as Jacobian
20 .  flag - matrix flag
21 
22    Options Database Key:
23 $  -snes_fd
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    SNESDefaultComputeJacobian() 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 .keywords: SNES, finite differences, Jacobian
33 
34 .seealso: SNESSetJacobian(), SNESTestJacobian()
35 @*/
36 int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,
37                                MatStructure *flag,void *ctx)
38 {
39   Vec    j1,j2,x2;
40   int    i,ierr,N,start,end,j;
41   Scalar dx, mone = -1.0,*y,scale,*xx;
42   double epsilon = 1.e-8,amax; /* assumes double precision */
43 
44   MatZeroEntries(*J);
45   ierr = VecDuplicate(x1,&j1); CHKERR(ierr);
46   ierr = VecDuplicate(x1,&j2); CHKERR(ierr);
47   ierr = VecDuplicate(x1,&x2); CHKERR(ierr);
48 
49   ierr = VecGetSize(x1,&N); CHKERR(ierr);
50   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERR(ierr);
51   VecGetArray(x1,&xx);
52   ierr = SNESComputeFunction(snes,x1,j1); CHKERR(ierr);
53   for ( i=0; i<N; i++ ) {
54     ierr = VecCopy(x1,x2); CHKERR(ierr);
55     if ( i>= start && i<end) {
56       dx = xx[i-start];
57       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
58       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
59       dx *= epsilon;
60       scale = -1.0/dx;
61       VecSetValues(x2,1,&i,&dx,ADDVALUES);
62     }
63     ierr = SNESComputeFunction(snes,x2,j2); CHKERR(ierr);
64     ierr = VecAXPY(&mone,j1,j2); CHKERR(ierr);
65     VecScale(&scale,j2);
66     VecGetArray(j2,&y);
67     VecAMax(j2,0,&amax); amax *= 1.e-14;
68     for ( j=start; j<end; j++ ) {
69       if (y[j-start] > amax || y[j-start] < -amax) {
70         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERR(ierr);
71       }
72     }
73     VecRestoreArray(j2,&y);
74   }
75   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
76   VecDestroy(x2); VecDestroy(j1); VecDestroy(j2);
77   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
78   return 0;
79 }
80 
81