xref: /petsc/src/snes/interface/snesj.c (revision 72ce74bd30d5d09d1371a3eb387704f4af395d1a)
1 
2 #include "src/snes/snesimpl.h"    /*I  "petscsnes.h"  I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESDefaultComputeJacobian"
6 /*@C
7    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
8 
9    Collective on SNES
10 
11    Input Parameters:
12 +  x1 - compute Jacobian at this point
13 -  ctx - application's function context, as set with SNESSetFunction()
14 
15    Output Parameters:
16 +  J - Jacobian matrix (not altered in this routine)
17 .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
18 -  flag - flag indicating whether the matrix sparsity structure has changed
19 
20    Options Database Key:
21 +  -snes_fd - Activates SNESDefaultComputeJacobian()
22 -  -snes_test_err - Square root of function error tolerance, default square root of machine
23                     epsilon (1.e-8 in double, 3.e-4 in single)
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    An alternative routine that uses coloring to explot matrix sparsity is
33    SNESDefaultComputeJacobianColor().
34 
35    Level: intermediate
36 
37 .keywords: SNES, finite differences, Jacobian
38 
39 .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor()
40 @*/
41 PetscErrorCode SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
42 {
43   Vec         j1a,j2a,x2;
44   int         i,ierr,N,start,end,j;
45   PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
46   PetscReal   amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
47   PetscReal   dx_min = 1.e-16,dx_par = 1.e-1;
48   MPI_Comm    comm;
49   int         (*eval_fct)(SNES,Vec,Vec)=0;
50 
51   PetscFunctionBegin;
52   ierr = PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr);
53   eval_fct = SNESComputeFunction;
54 
55   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
56   ierr = MatZeroEntries(*B);CHKERRQ(ierr);
57   if (!snes->nvwork) {
58     ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr);
59     snes->nvwork = 3;
60     PetscLogObjectParents(snes,3,snes->vwork);
61   }
62   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
63 
64   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
65   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
66   ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr);
67 
68   /* Compute Jacobian approximation, 1 column at a time.
69       x1 = current iterate, j1a = F(x1)
70       x2 = perturbed iterate, j2a = F(x2)
71    */
72   for (i=0; i<N; i++) {
73     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
74     if (i>= start && i<end) {
75       ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
76       dx = xx[i-start];
77       ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
78 #if !defined(PETSC_USE_COMPLEX)
79       if (dx < dx_min && dx >= 0.0) dx = dx_par;
80       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
81 #else
82       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
83       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
84 #endif
85       dx *= epsilon;
86       wscale = 1.0/dx;
87       ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
88     } else {
89       wscale = 0.0;
90     }
91     ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr);
92     ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr);
93     /* Communicate scale to all processors */
94     ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
95     ierr = VecScale(&scale,j2a);CHKERRQ(ierr);
96     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
97     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
98     for (j=start; j<end; j++) {
99       if (PetscAbsScalar(y[j-start]) > amax) {
100         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
101       }
102     }
103     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
104   }
105   ierr  = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106   ierr  = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107   *flag =  DIFFERENT_NONZERO_PATTERN;
108   PetscFunctionReturn(0);
109 }
110 
111 
112