xref: /petsc/src/snes/interface/snesj.c (revision da88328c9011629587562b3a985d9f317f84061b)
1 
2 
3 #ifndef lint
4 static char vcid[] = "$Id: snesj.c,v 1.10 1995/05/12 04:18:43 bsmith Exp curfman $";
5 #endif
6 
7 #include "draw.h"
8 #include "snes.h"
9 #include "options.h"
10 
11 /*@
12    SNESDefaultComputeJacobian - Computes the Jacobian using finite
13    differences.
14 
15    Input Parameters:
16 .  x1 - compute Jacobian at this point
17 .  ctx - application's function context, as set with SNESSetFunction()
18 
19    Output Parameters:
20 .  J - Jacobian
21 .  B - preconditioner, same as Jacobian
22 .  flag - matrix flag
23 
24    Notes:
25    This routine is slow and expensive, and is not currently optimized
26    to take advantage of sparsity in the problem.  Although
27    SNESDefaultComputeJacobian() is not recommended for general use
28    in large-scale applications, It can be useful in checking the
29    correctness of a user-provided Jacobian.
30 
31 .keywords: SNES, finite differences, Jacobian
32 
33 .seealso: SNESSetJacobian(), SNESTestJacobian()
34 @*/
35 int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,int *flag,
36                                void *ctx)
37 {
38   Vec    j1,j2,x2;
39   int    i,ierr,N,start,end,j;
40   Scalar dx, mone = -1.0,*y,scale,*xx;
41   double epsilon = 1.e-8,amax; /* assumes double precision */
42 
43   MatZeroEntries(*J);
44   ierr = VecDuplicate(x1,&j1); CHKERR(ierr);
45   ierr = VecDuplicate(x1,&j2); CHKERR(ierr);
46   ierr = VecDuplicate(x1,&x2); CHKERR(ierr);
47 
48   ierr = VecGetSize(x1,&N); CHKERR(ierr);
49   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERR(ierr);
50   VecGetArray(x1,&xx);
51   ierr = SNESComputeFunction(snes,x1,j1); CHKERR(ierr);
52   for ( i=0; i<N; i++ ) {
53     ierr = VecCopy(x1,x2); CHKERR(ierr);
54     if ( i>= start && i<end) {
55       dx = xx[i-start];
56       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
57       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
58       dx *= epsilon;
59       scale = -1.0/dx;
60       VecSetValues(x2,1,&i,&dx,ADDVALUES);
61     }
62     ierr = SNESComputeFunction(snes,x2,j2); CHKERR(ierr);
63     ierr = VecAXPY(&mone,j1,j2); CHKERR(ierr);
64     VecScale(&scale,j2);
65     VecGetArray(j2,&y);
66     VecAMax(j2,0,&amax); amax *= 1.e-14;
67     for ( j=start; j<end; j++ ) {
68       if (y[j-start] > amax || y[j-start] < -amax) {
69         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERR(ierr);
70       }
71     }
72     VecRestoreArray(j2,&y);
73   }
74   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
75   VecDestroy(x2); VecDestroy(j1); VecDestroy(j2);
76   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
77   return 0;
78 }
79 
80