Actual source code: snesj.c

  1: /*$Id: snesj.c,v 1.75 2001/09/11 18:06:40 bsmith Exp $*/

 3:  #include src/snes/snesimpl.h

  5: /*@C
  6:    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 

  8:    Collective on SNES

 10:    Input Parameters:
 11: +  x1 - compute Jacobian at this point
 12: -  ctx - application's function context, as set with SNESSetFunction()

 14:    Output Parameters:
 15: +  J - Jacobian matrix (not altered in this routine)
 16: .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
 17: -  flag - flag indicating whether the matrix sparsity structure has changed

 19:    Options Database Key:
 20: +  -snes_fd - Activates SNESDefaultComputeJacobian()
 21: -  -snes_test_err - Square root of function error tolerance, default square root of machine
 22:                     epsilon (1.e-8 in double, 3.e-4 in single)

 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.

 31:    An alternative routine that uses coloring to explot matrix sparsity is
 32:    SNESDefaultComputeJacobianColor().

 34:    Level: intermediate

 36: .keywords: SNES, finite differences, Jacobian

 38: .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor()
 39: @*/
 40: int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
 41: {
 42:   Vec         j1a,j2a,x2;
 43:   int         i,ierr,N,start,end,j;
 44:   PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
 45:   PetscReal   amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
 46:   PetscReal   dx_min = 1.e-16,dx_par = 1.e-1;
 47:   MPI_Comm    comm;
 48:   int         (*eval_fct)(SNES,Vec,Vec)=0;

 51:   PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);
 52:   eval_fct = SNESComputeFunction;

 54:   PetscObjectGetComm((PetscObject)x1,&comm);
 55:   MatZeroEntries(*B);
 56:   if (!snes->nvwork) {
 57:     VecDuplicateVecs(x1,3,&snes->vwork);
 58:     snes->nvwork = 3;
 59:     PetscLogObjectParents(snes,3,snes->vwork);
 60:   }
 61:   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];

 63:   VecGetSize(x1,&N);
 64:   VecGetOwnershipRange(x1,&start,&end);
 65:   (*eval_fct)(snes,x1,j1a);

 67:   /* Compute Jacobian approximation, 1 column at a time. 
 68:       x1 = current iterate, j1a = F(x1)
 69:       x2 = perturbed iterate, j2a = F(x2)
 70:    */
 71:   for (i=0; i<N; i++) {
 72:     VecCopy(x1,x2);
 73:     if (i>= start && i<end) {
 74:       VecGetArray(x1,&xx);
 75:       dx = xx[i-start];
 76:       VecRestoreArray(x1,&xx);
 77: #if !defined(PETSC_USE_COMPLEX)
 78:       if (dx < dx_min && dx >= 0.0) dx = dx_par;
 79:       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
 80: #else
 81:       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
 82:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
 83: #endif
 84:       dx *= epsilon;
 85:       wscale = 1.0/dx;
 86:       VecSetValues(x2,1,&i,&dx,ADD_VALUES);
 87:     } else {
 88:       wscale = 0.0;
 89:     }
 90:     (*eval_fct)(snes,x2,j2a);
 91:     VecAXPY(&mone,j1a,j2a);
 92:     /* Communicate scale to all processors */
 93:     MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);
 94:     VecScale(&scale,j2a);
 95:     VecNorm(j2a,NORM_INFINITY,&amax); amax *= 1.e-14;
 96:     VecGetArray(j2a,&y);
 97:     for (j=start; j<end; j++) {
 98:       if (PetscAbsScalar(y[j-start]) > amax) {
 99:         MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);
100:       }
101:     }
102:     VecRestoreArray(j2a,&y);
103:   }
104:   ierr  = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
105:   ierr  = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
106:   *flag =  DIFFERENT_NONZERO_PATTERN;
107:   return(0);
108: }