xref: /petsc/src/snes/interface/snesj.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
111320018SBarry Smith 
2e090d566SSatish Balay #include "src/snes/snesimpl.h"    /*I  "petscsnes.h"  I*/
311320018SBarry Smith 
44a2ae208SSatish Balay #undef __FUNCT__
54a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultComputeJacobian"
64b828684SBarry Smith /*@C
784a7bf8dSLois Curfman McInnes    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
811320018SBarry Smith 
9fee21e36SBarry Smith    Collective on SNES
10fee21e36SBarry Smith 
11c7afd0dbSLois Curfman McInnes    Input Parameters:
12c7afd0dbSLois Curfman McInnes +  x1 - compute Jacobian at this point
13c7afd0dbSLois Curfman McInnes -  ctx - application's function context, as set with SNESSetFunction()
14c7afd0dbSLois Curfman McInnes 
15c7afd0dbSLois Curfman McInnes    Output Parameters:
16c7afd0dbSLois Curfman McInnes +  J - Jacobian matrix (not altered in this routine)
17c7afd0dbSLois Curfman McInnes .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
18c7afd0dbSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
19c7afd0dbSLois Curfman McInnes 
20ad960d00SLois Curfman McInnes    Options Database Key:
2108f75eefSBarry Smith +  -snes_fd - Activates SNESDefaultComputeJacobian()
2277d8c4bbSBarry Smith -  -snes_test_err - Square root of function error tolerance, default square root of machine
2377d8c4bbSBarry Smith                     epsilon (1.e-8 in double, 3.e-4 in single)
24ad960d00SLois Curfman McInnes 
255f3c43d9SLois Curfman McInnes    Notes:
265f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
275f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
285f3c43d9SLois Curfman McInnes    SNESDefaultComputeJacobian() is not recommended for general use
295f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
305f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3111320018SBarry Smith 
32b4fc646aSLois Curfman McInnes    An alternative routine that uses coloring to explot matrix sparsity is
332d0c0e3bSBarry Smith    SNESDefaultComputeJacobianColor().
34b4fc646aSLois Curfman McInnes 
3536851e7fSLois Curfman McInnes    Level: intermediate
3636851e7fSLois Curfman McInnes 
375f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian
385f3c43d9SLois Curfman McInnes 
392d0c0e3bSBarry Smith .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor()
4011320018SBarry Smith @*/
41dfbe8321SBarry Smith PetscErrorCode SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
4211320018SBarry Smith {
4388c956adSLois Curfman McInnes   Vec         j1a,j2a,x2;
44*6849ba73SBarry Smith   PetscErrorCode ierr;
45*6849ba73SBarry Smith   int         i,N,start,end,j;
46ea709b57SSatish Balay   PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
4777d8c4bbSBarry Smith   PetscReal   amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
48329f5518SBarry Smith   PetscReal   dx_min = 1.e-16,dx_par = 1.e-1;
49bbb6d6a8SBarry Smith   MPI_Comm    comm;
50*6849ba73SBarry Smith   PetscErrorCode (*eval_fct)(SNES,Vec,Vec)=0;
510521c3abSLois Curfman McInnes 
523a40ed3dSBarry Smith   PetscFunctionBegin;
5387828ca2SBarry Smith   ierr = PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr);
544b27c08aSLois Curfman McInnes   eval_fct = SNESComputeFunction;
5523242f5aSBarry Smith 
562d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
572d0c0e3bSBarry Smith   ierr = MatZeroEntries(*B);CHKERRQ(ierr);
58aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
59aa79bc6dSLois Curfman McInnes     ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr);
60aa79bc6dSLois Curfman McInnes     snes->nvwork = 3;
61b0a32e0cSBarry Smith     PetscLogObjectParents(snes,3,snes->vwork);
62aa79bc6dSLois Curfman McInnes   }
6388c956adSLois Curfman McInnes   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
6423242f5aSBarry Smith 
6578b31e54SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
6678b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
67a86fb38aSBarry Smith   ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr);
68c005e166SLois Curfman McInnes 
69c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
7088c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
7188c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
72c005e166SLois Curfman McInnes    */
7339e2f89bSBarry Smith   for (i=0; i<N; i++) {
7478b31e54SBarry Smith     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
7523242f5aSBarry Smith     if (i>= start && i<end) {
762e8a6d31SBarry Smith       ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
7739e2f89bSBarry Smith       dx = xx[i-start];
782e8a6d31SBarry Smith       ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
79aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8054d73c9fSLois Curfman McInnes       if (dx < dx_min && dx >= 0.0) dx = dx_par;
8154d73c9fSLois Curfman McInnes       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
8219a167f6SBarry Smith #else
83329f5518SBarry Smith       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
84329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
8519a167f6SBarry Smith #endif
8639e2f89bSBarry Smith       dx *= epsilon;
8774f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
882e8a6d31SBarry Smith       ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
896c783aadSBarry Smith     } else {
90bbb6d6a8SBarry Smith       wscale = 0.0;
91bbb6d6a8SBarry Smith     }
92a86fb38aSBarry Smith     ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr);
9388c956adSLois Curfman McInnes     ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr);
94c005e166SLois Curfman McInnes     /* Communicate scale to all processors */
956831982aSBarry Smith     ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
962e8a6d31SBarry Smith     ierr = VecScale(&scale,j2a);CHKERRQ(ierr);
972e8a6d31SBarry Smith     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
98d142edc7SBarry Smith     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
9923242f5aSBarry Smith     for (j=start; j<end; j++) {
100cddf8d76SBarry Smith       if (PetscAbsScalar(y[j-start]) > amax) {
101b4fc646aSLois Curfman McInnes         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
10223242f5aSBarry Smith       }
10323242f5aSBarry Smith     }
1042e8a6d31SBarry Smith     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
10523242f5aSBarry Smith   }
106b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108595b82d2SBarry Smith   *flag =  DIFFERENT_NONZERO_PATTERN;
1093a40ed3dSBarry Smith   PetscFunctionReturn(0);
11011320018SBarry Smith }
11111320018SBarry Smith 
112bd45824fSLois Curfman McInnes 
113