Actual source code: snesj.c

  1: #define PETSCSNES_DLL

 3:  #include include/private/snesimpl.h

  7: /*@C
  8:    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 

 10:    Collective on SNES

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

 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

 21:    Options Database Key:
 22: +  -snes_fd - Activates SNESDefaultComputeJacobian()
 23: .  -snes_test_err - Square root of function error tolerance, default square root of machine
 24:                     epsilon (1.e-8 in double, 3.e-4 in single)
 25: -  -mat_fd_type - Either wp or ds (see MATSNESMF_WP or MATSNESMF_DS)

 27:    Notes:
 28:    This routine is slow and expensive, and is not currently optimized
 29:    to take advantage of sparsity in the problem.  Although
 30:    SNESDefaultComputeJacobian() is not recommended for general use
 31:    in large-scale applications, It can be useful in checking the
 32:    correctness of a user-provided Jacobian.

 34:    An alternative routine that uses coloring to exploit matrix sparsity is
 35:    SNESDefaultComputeJacobianColor().

 37:    Level: intermediate

 39: .keywords: SNES, finite differences, Jacobian

 41: .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor(), MatCreateSNESMF()
 42: @*/
 43: PetscErrorCode  SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
 44: {
 45:   Vec            j1a,j2a,x2;
 47:   PetscInt       i,N,start,end,j,value;
 48:   PetscScalar    dx,*y,scale,*xx,wscale;
 49:   PetscReal      amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
 50:   PetscReal      dx_min = 1.e-16,dx_par = 1.e-1,unorm;
 51:   MPI_Comm       comm;
 52:   PetscErrorCode (*eval_fct)(SNES,Vec,Vec)=0;
 53:   PetscTruth     assembled,use_wp = PETSC_TRUE,flg;
 54:   const char     *list[2] = {"ds","wp"};

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

 60:   PetscObjectGetComm((PetscObject)x1,&comm);
 61:   MatAssembled(*B,&assembled);
 62:   if (assembled) {
 63:     MatZeroEntries(*B);
 64:   }
 65:   if (!snes->nvwork) {
 66:     VecDuplicateVecs(x1,3,&snes->vwork);
 67:     snes->nvwork = 3;
 68:     PetscLogObjectParents(snes,3,snes->vwork);
 69:   }
 70:   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];

 72:   VecGetSize(x1,&N);
 73:   VecGetOwnershipRange(x1,&start,&end);
 74:   (*eval_fct)(snes,x1,j1a);

 76:   PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESDefaultComputeJacobian",list,2,"wp",&value,&flg);
 77:   if (flg && !value) {
 78:     use_wp = PETSC_FALSE;
 79:   }
 80:   if (use_wp) {
 81:     VecNorm(x1,NORM_2,&unorm);
 82:   }
 83:   /* Compute Jacobian approximation, 1 column at a time. 
 84:       x1 = current iterate, j1a = F(x1)
 85:       x2 = perturbed iterate, j2a = F(x2)
 86:    */
 87:   for (i=0; i<N; i++) {
 88:     VecCopy(x1,x2);
 89:     if (i>= start && i<end) {
 90:       VecGetArray(x1,&xx);
 91:       if (use_wp) {
 92:         dx = 1.0 + unorm;
 93:       } else {
 94:         dx = xx[i-start];
 95:       }
 96:       VecRestoreArray(x1,&xx);
 97: #if !defined(PETSC_USE_COMPLEX)
 98:       if (dx < dx_min && dx >= 0.0) dx = dx_par;
 99:       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
100: #else
101:       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
102:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
103: #endif
104:       dx *= epsilon;
105:       wscale = 1.0/dx;
106:       VecSetValues(x2,1,&i,&dx,ADD_VALUES);
107:     } else {
108:       wscale = 0.0;
109:     }
110:     (*eval_fct)(snes,x2,j2a);
111:     VecAXPY(j2a,-1.0,j1a);
112:     /* Communicate scale to all processors */
113:     MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);
114:     VecScale(j2a,scale);
115:     VecNorm(j2a,NORM_INFINITY,&amax); amax *= 1.e-14;
116:     VecGetArray(j2a,&y);
117:     for (j=start; j<end; j++) {
118:       if (PetscAbsScalar(y[j-start]) > amax) {
119:         MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);
120:       }
121:     }
122:     VecRestoreArray(j2a,&y);
123:   }
124:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
125:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
126:   if (*B != *J) {
127:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
128:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
129:   }
130:   *flag =  DIFFERENT_NONZERO_PATTERN;
131:   return(0);
132: }