xref: /petsc/src/snes/interface/snesj.c (revision 6831982abb6453c2f3e25efecb5d0951d333371e)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: snesj.c,v 1.60 1999/10/22 21:16:38 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    Collective on SNES
13 
14    Input Parameters:
15 +  x1 - compute Jacobian at this point
16 -  ctx - application's function context, as set with SNESSetFunction()
17 
18    Output Parameters:
19 +  J - Jacobian matrix (not altered in this routine)
20 .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
21 -  flag - flag indicating whether the matrix sparsity structure has changed
22 
23    Options Database Key:
24 .  -snes_fd - Activates SNESDefaultComputeJacobian()
25 
26    Notes:
27    This routine is slow and expensive, and is not currently optimized
28    to take advantage of sparsity in the problem.  Although
29    SNESDefaultComputeJacobian() is not recommended for general use
30    in large-scale applications, It can be useful in checking the
31    correctness of a user-provided Jacobian.
32 
33    An alternative routine that uses coloring to explot matrix sparsity is
34    SNESDefaultComputeJacobianColor().
35 
36    Level: intermediate
37 
38 .keywords: SNES, finite differences, Jacobian
39 
40 .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor()
41 @*/
42 int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
43 {
44   Vec      j1a,j2a,x2;
45   int      i,ierr,N,start,end,j;
46   Scalar   dx, mone = -1.0,*y,scale,*xx,wscale;
47   double   amax, epsilon = 1.e-8; /* assumes double precision */
48   double   dx_min = 1.e-16, dx_par = 1.e-1;
49   MPI_Comm comm;
50   int      (*eval_fct)(SNES,Vec,Vec)=0;
51 
52   PetscFunctionBegin;
53   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction;
54   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient;
55   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
56 
57   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
58   ierr = MatZeroEntries(*B);CHKERRQ(ierr);
59   if (!snes->nvwork) {
60     ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr);
61     snes->nvwork = 3;
62     PLogObjectParents(snes,3,snes->vwork);
63   }
64   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
65 
66   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
67   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
68   ierr = eval_fct(snes,x1,j1a);CHKERRQ(ierr);
69 
70   /* Compute Jacobian approximation, 1 column at a time.
71       x1 = current iterate, j1a = F(x1)
72       x2 = perturbed iterate, j2a = F(x2)
73    */
74   for ( i=0; i<N; i++ ) {
75     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
76     if ( i>= start && i<end) {
77       ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
78       dx = xx[i-start];
79       ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
80 #if !defined(PETSC_USE_COMPLEX)
81       if (dx < dx_min && dx >= 0.0) dx = dx_par;
82       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
83 #else
84       if (PetscAbsScalar(dx) < dx_min && PetscReal(dx) >= 0.0) dx = dx_par;
85       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
86 #endif
87       dx *= epsilon;
88       wscale = 1.0/dx;
89       ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
90     } else {
91       wscale = 0.0;
92     }
93     ierr = eval_fct(snes,x2,j2a);CHKERRQ(ierr);
94     ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr);
95     /* Communicate scale to all processors */
96     ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
97     ierr = VecScale(&scale,j2a);CHKERRQ(ierr);
98     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
99     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); 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     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
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    Collective on SNES
119 
120    Input Parameters:
121 +  x1 - compute Hessian at this point
122 -  ctx - application's gradient context, as set with SNESSetGradient()
123 
124    Output Parameters:
125 +  J - Hessian matrix (not altered in this routine)
126 .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
127 -  flag - flag indicating whether the matrix sparsity structure has changed
128 
129    Options Database Key:
130 $  -snes_fd - Activates SNESDefaultComputeHessian()
131 
132 
133    Level: intermediate
134 
135    Notes:
136    This routine is slow and expensive, and is not currently optimized
137    to take advantage of sparsity in the problem.  Although
138    SNESDefaultComputeHessian() is not recommended for general use
139    in large-scale applications, It can be useful in checking the
140    correctness of a user-provided Hessian.
141 
142 .keywords: SNES, finite differences, Hessian
143 
144 .seealso: SNESSetHessian()
145 @*/
146 int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,
147                               MatStructure *flag,void *ctx)
148 {
149   int ierr;
150 
151   PetscFunctionBegin;
152   ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNC__
157 #define __FUNC__ "SNESDefaultComputeHessianColor"
158 /*@C
159    SNESDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
160 
161    Collective on SNES
162 
163    Input Parameters:
164 +  x1 - compute Hessian at this point
165 -  ctx - application's gradient context, as set with SNESSetGradient()
166 
167    Output Parameters:
168 +  J - Hessian matrix (not altered in this routine)
169 .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
170 -  flag - flag indicating whether the matrix sparsity structure has changed
171 
172     Options Database Keys:
173 .  -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor()
174 
175    Level: intermediate
176 
177  .keywords: SNES, finite differences, Hessian, coloring, sparse
178 
179 .seealso: SNESSetHessian()
180 @*/
181 int SNESDefaultComputeHessianColor(SNES snes,Vec x1,Mat *J,Mat *B,
182                               MatStructure *flag,void *ctx)
183 {
184   int ierr;
185 
186   PetscFunctionBegin;
187   ierr = SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
188   PetscFunctionReturn(0);
189 }
190 
191