xref: /petsc/src/snes/interface/snesj.c (revision aa79bc6d89adb718802453cbc2f155d610aa3690)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: snesj.c,v 1.27 1996/02/15 02:06:09 curfman Exp curfman $";
4 #endif
5 
6 #include "draw.h"    /*I  "draw.h"  I*/
7 #include "snesimpl.h"    /*I  "snes.h"  I*/
8 
9 /*@C
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,MatStructure *flag,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,wscale;
41   double   amax, epsilon = 1.e-8; /* assumes double precision */
42   MPI_Comm comm;
43   int      (*eval_fct)(SNES,Vec,Vec);
44 
45   if (snes->method_class == SNES_NONLINEAR_EQUATIONS)
46     eval_fct = SNESComputeFunction;
47   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
48     eval_fct = SNESComputeGradient;
49   else SETERRQ(1,"SNESDefaultComputeJacobian: Invalid method class");
50 
51   PetscObjectGetComm((PetscObject)x1,&comm);
52   MatZeroEntries(*J);
53   if (!snes->nvwork) {
54     ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr);
55     snes->nvwork = 3;
56     PLogObjectParents(snes,3,snes->vwork);
57   }
58   j1 = snes->vwork[0]; j2 = snes->vwork[1]; x2 = snes->vwork[2];
59 
60   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
61   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
62   VecGetArray(x1,&xx);
63   ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr);
64   for ( i=0; i<N; i++ ) {
65     ierr = VecCopy(x1,x2); CHKERRQ(ierr);
66     if ( i>= start && i<end) {
67       dx = xx[i-start];
68 #if !defined(PETSC_COMPLEX)
69       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
70       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
71 #else
72       if (abs(dx) < 1.e-16 && real(dx) >= 0.0) dx = 1.e-1;
73       else if (real(dx) < 0.0 && abs(dx) < 1.e-16) dx = -1.e-1;
74 #endif
75       dx *= epsilon;
76       wscale = -1.0/dx;
77       VecSetValues(x2,1,&i,&dx,ADD_VALUES);
78     }
79     else {
80       wscale = 0.0;
81     }
82     ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr);
83     ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr);
84 /* communicate scale to all processors */
85 #if !defined(PETSC_COMPLEX)
86     MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);
87 #else
88     MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);
89 #endif
90     scale = -scale;
91     VecScale(&scale,j2);
92     VecGetArray(j2,&y);
93     VecNorm(j2,NORM_INFINITY,&amax); amax *= 1.e-14;
94     for ( j=start; j<end; j++ ) {
95       if (PetscAbsScalar(y[j-start]) > amax) {
96         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr);
97       }
98     }
99     VecRestoreArray(j2,&y);
100   }
101   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
102   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
103   return 0;
104 }
105 
106 /*@C
107    SNESDefaultComputeHessian - Computes the Hessian using finite
108    differences.
109 
110    Input Parameters:
111 .  x1 - compute Hessian at this point
112 .  ctx - application's gradient context, as set with SNESSetGradient()
113 
114    Output Parameters:
115 .  J - Hessian
116 .  B - preconditioner, same as Hessian
117 .  flag - matrix flag
118 
119    Options Database Key:
120 $  -snes_fd
121 
122    Notes:
123    This routine is slow and expensive, and is not currently optimized
124    to take advantage of sparsity in the problem.  Although
125    SNESDefaultComputeHessian() is not recommended for general use
126    in large-scale applications, It can be useful in checking the
127    correctness of a user-provided Hessian.
128 
129 .keywords: SNES, finite differences, Hessian
130 
131 .seealso: SNESSetHessian(), SNESTestHessian()
132 @*/
133 int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
134 {
135   return SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);
136 }
137