xref: /petsc/src/snes/interface/snesj.c (revision 336a5e98162d189faf560c36b5de145b11386e01)
1 
2 
3 #ifndef lint
4 static char vcid[] = "$Id: snesj.c,v 1.7 1995/05/05 03:51:29 bsmith Exp bsmith $";
5 #endif
6 
7 #include "draw.h"
8 #include "snes.h"
9 #include "options.h"
10 
11 /*@
12    SNESDefaultComputeJacobian - Computes Jacobian using finite
13        differences. Slow and expensive.
14 
15  Input Parameters:
16 .  x - compute Jacobian at this point
17 .  ctx - applications Function context
18 
19   Output Parameters:
20 .  J - Jacobian
21 .  B - preconditioner, same as Jacobian
22 
23 .keywords: finite differences, Jacobian
24 
25 .seealso: SNESSetJacobian, SNESTestJacobian
26 @*/
27 int SNESDefaultComputeJacobian(SNES snes, Vec x1,Mat *J,Mat *B,int *flag,
28                                void *ctx)
29 {
30   Vec    j1,j2,x2;
31   int    i,ierr,N,start,end,j;
32   Scalar dx, mone = -1.0,*y,scale,*xx;
33   double epsilon = 1.e-8,amax; /* assumes double precision */
34 
35   MatZeroEntries(*J);
36   ierr = VecDuplicate(x1,&j1); CHKERR(ierr);
37   ierr = VecDuplicate(x1,&j2); CHKERR(ierr);
38   ierr = VecDuplicate(x1,&x2); CHKERR(ierr);
39 
40   ierr = VecGetSize(x1,&N); CHKERR(ierr);
41   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERR(ierr);
42   VecGetArray(x1,&xx);
43   ierr = SNESComputeFunction(snes,x1,j1); CHKERR(ierr);
44   for ( i=0; i<N; i++ ) {
45     ierr = VecCopy(x1,x2); CHKERR(ierr);
46     if ( i>= start && i<end) {
47       dx = xx[i-start];
48       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
49       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
50       dx *= epsilon;
51       scale = -1.0/dx;
52       VecSetValues(x2,1,&i,&dx,ADDVALUES);
53     }
54     ierr = SNESComputeFunction(snes,x2,j2); CHKERR(ierr);
55     ierr = VecAXPY(&mone,j1,j2); CHKERR(ierr);
56     VecScale(&scale,j2);
57     VecGetArray(j2,&y);
58     VecAMax(j2,0,&amax); amax *= 1.e-14;
59     for ( j=start; j<end; j++ ) {
60       if (y[j-start] > amax || y[j-start] < -amax) {
61         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERR(ierr);
62       }
63     }
64     VecRestoreArray(j2,&y);
65   }
66   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
67   VecDestroy(x2); VecDestroy(j1); VecDestroy(j2);
68   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
69   return 0;
70 }
71 
72