xref: /petsc/src/snes/interface/snesj.c (revision c01c455d43b147fbcccdf73d28be4d01ad7032e0)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: snesj.c,v 1.19 1995/08/14 19:39:37 bsmith Exp curfman $";
4 #endif
5 
6 #include "draw.h"    /*I  "draw.h"  I*/
7 #include "snes.h"    /*I  "snes.h"  I*/
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,wscale;
42   double   epsilon = 1.e-8,amax; /* assumes double precision */
43   MPI_Comm comm;
44 
45   PetscObjectGetComm((PetscObject)x1,&comm);
46   MatZeroEntries(*J);
47   ierr = VecDuplicate(x1,&j1); CHKERRQ(ierr);
48   ierr = VecDuplicate(x1,&j2); CHKERRQ(ierr);
49   ierr = VecDuplicate(x1,&x2); CHKERRQ(ierr);
50   PLogObjectParent(snes,j1); PLogObjectParent(snes,j2);
51   PLogObjectParent(snes,x2);
52 
53   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
54   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
55   VecGetArray(x1,&xx);
56   ierr = SNESComputeFunction(snes,x1,j1); CHKERRQ(ierr);
57   for ( i=0; i<N; i++ ) {
58     ierr = VecCopy(x1,x2); CHKERRQ(ierr);
59     if ( i>= start && i<end) {
60       dx = xx[i-start];
61 #if !defined(PETSC_COMPLEX)
62       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
63       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
64 #else
65       if (abs(dx) < 1.e-16 && real(dx) >= 0.0) dx = 1.e-1;
66       else if (real(dx) < 0.0 && abs(dx) > -1.e-16) dx = -1.e-1;
67 #endif
68       dx *= epsilon;
69       wscale = -1.0/dx;
70       VecSetValues(x2,1,&i,&dx,ADDVALUES);
71     }
72     else {
73       wscale = 0.0;
74     }
75     ierr = SNESComputeFunction(snes,x2,j2); CHKERRQ(ierr);
76     ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr);
77 /* communicate scale to all processors */
78 #if !defined(PETSC_COMPLEX)
79     MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);
80 #else
81     MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);
82 #endif
83     VecScale(&scale,j2);
84     VecGetArray(j2,&y);
85     VecAMax(j2,0,&amax); amax *= 1.e-14;
86     for ( j=start; j<end; j++ ) {
87 #if defined(PETSC_COMPLEX)
88       if (abs(y[j-start]) > amax) {
89 #else
90       if (y[j-start] > amax || y[j-start] < -amax) {
91 #endif
92         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERRQ(ierr);
93       }
94     }
95     VecRestoreArray(j2,&y);
96   }
97   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
98   VecDestroy(x2); VecDestroy(j1); VecDestroy(j2);
99   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
100   return 0;
101 }
102 
103