xref: /petsc/src/snes/interface/snesj.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: snesj.c,v 1.44 1997/09/25 22:42:32 curfman Exp bsmith $";
3 #endif
4 
5 #include "src/snes/snesimpl.h"    /*I  "snes.h"  I*/
6 
7 #undef __FUNC__
8 #define __FUNC__ "SNESDefaultComputeJacobian"
9 /*@C
10    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
11 
12    Input Parameters:
13 .  x1 - compute Jacobian at this point
14 .  ctx - application's function context, as set with SNESSetFunction()
15 
16    Output Parameters:
17 .  J - Jacobian matrix (not altered in this routine)
18 .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
19 .  flag - flag indicating whether the matrix sparsity structure has changed
20 
21    Options Database Key:
22 $  -snes_fd
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    An alternative routine that uses coloring to explot matrix sparsity is
32    SNESDefaultComputeJacobianWithColoring().
33 
34 .keywords: SNES, finite differences, Jacobian
35 
36 .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianWithColoring()
37 @*/
38 int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
39 {
40   Vec      j1,j2,x2;
41   int      i,ierr,N,start,end,j;
42   Scalar   dx, mone = -1.0,*y,scale,*xx,wscale;
43   double   amax, epsilon = 1.e-8; /* assumes double precision */
44   double   dx_min = 1.e-16, dx_par = 1.e-1;
45   MPI_Comm comm;
46   int      (*eval_fct)(SNES,Vec,Vec);
47 
48   PetscFunctionBegin;
49   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction;
50   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient;
51   else SETERRQ(1,0,"Invalid method class");
52 
53   PetscObjectGetComm((PetscObject)x1,&comm);
54   MatZeroEntries(*B);
55   if (!snes->nvwork) {
56     ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr);
57     snes->nvwork = 3;
58     PLogObjectParents(snes,3,snes->vwork);
59   }
60   j1 = snes->vwork[0]; j2 = snes->vwork[1]; x2 = snes->vwork[2];
61 
62   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
63   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
64   VecGetArray(x1,&xx);
65   ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr);
66 
67   /* Compute Jacobian approximation, 1 column at a time.
68       x1 = current iterate, j1 = F(x1)
69       x2 = perturbed iterate, j2 = F(x2)
70    */
71   for ( i=0; i<N; i++ ) {
72     ierr = VecCopy(x1,x2); CHKERRQ(ierr);
73     if ( i>= start && i<end) {
74       dx = xx[i-start];
75 #if !defined(USE_PETSC_COMPLEX)
76       if (dx < dx_min && dx >= 0.0) dx = dx_par;
77       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
78 #else
79       if (abs(dx) < dx_min && real(dx) >= 0.0) dx = dx_par;
80       else if (real(dx) < 0.0 && abs(dx) < dx_min) dx = -dx_par;
81 #endif
82       dx *= epsilon;
83       wscale = 1.0/dx;
84       VecSetValues(x2,1,&i,&dx,ADD_VALUES);
85     }
86     else {
87       wscale = 0.0;
88     }
89     ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr);
90     ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr);
91     /* Communicate scale to all processors */
92 #if !defined(USE_PETSC_COMPLEX)
93     MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);
94 #else
95     MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);
96 #endif
97     VecScale(&scale,j2);
98     VecGetArray(j2,&y);
99     VecNorm(j2,NORM_INFINITY,&amax); amax *= 1.e-14;
100     for ( j=start; j<end; j++ ) {
101       if (PetscAbsScalar(y[j-start]) > amax) {
102         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr);
103       }
104     }
105     VecRestoreArray(j2,&y);
106   }
107   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
108   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
109   *flag =  DIFFERENT_NONZERO_PATTERN;
110   PetscFunctionReturn(0);
111 }
112 
113 #undef __FUNC__
114 #define __FUNC__ "SNESDefaultComputeHessian"
115 /*@C
116    SNESDefaultComputeHessian - Computes the Hessian using finite differences.
117 
118    Input Parameters:
119 .  x1 - compute Hessian at this point
120 .  ctx - application's gradient context, as set with SNESSetGradient()
121 
122    Output Parameters:
123 .  J - Hessian matrix (not altered in this routine)
124 .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
125 .  flag - flag indicating whether the matrix sparsity structure has changed
126 
127    Options Database Key:
128 $  -snes_fd
129 
130    Notes:
131    This routine is slow and expensive, and is not currently optimized
132    to take advantage of sparsity in the problem.  Although
133    SNESDefaultComputeHessian() is not recommended for general use
134    in large-scale applications, It can be useful in checking the
135    correctness of a user-provided Hessian.
136 
137 .keywords: SNES, finite differences, Hessian
138 
139 .seealso: SNESSetHessian()
140 @*/
141 int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
142 {
143   int ierr;
144 
145   PetscFunctionBegin;
146   ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149