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:   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction;
 53:   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient;
 54:   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid method class");

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

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

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

112: /*@C
113:    SNESDefaultComputeHessian - Computes the Hessian using finite differences. 

115:    Collective on SNES

117:    Input Parameters:
118: +  x1 - compute Hessian at this point
119: -  ctx - application's gradient context, as set with SNESSetGradient()

121:    Output Parameters:
122: +  J - Hessian matrix (not altered in this routine)
123: .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
124: -  flag - flag indicating whether the matrix sparsity structure has changed

126:    Options Database Key:
127: $  -snes_fd - Activates SNESDefaultComputeHessian()


130:    Level: intermediate

132:    Notes:
133:    This routine is slow and expensive, and is not currently optimized
134:    to take advantage of sparsity in the problem.  Although
135:    SNESDefaultComputeHessian() is not recommended for general use
136:    in large-scale applications, It can be useful in checking the
137:    correctness of a user-provided Hessian.

139: .keywords: SNES, finite differences, Hessian

141: .seealso: SNESSetHessian()
142: @*/
143: int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
144: {

148:   SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);
149:   return(0);
150: }

152: /*@C
153:    SNESDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 

155:    Collective on SNES

157:    Input Parameters:
158: +  x1 - compute Hessian at this point
159: -  ctx - application's gradient context, as set with SNESSetGradient()

161:    Output Parameters:
162: +  J - Hessian matrix (not altered in this routine)
163: .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
164: -  flag - flag indicating whether the matrix sparsity structure has changed

166:     Options Database Keys:
167: .  -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor()

169:    Level: intermediate

171:  .keywords: SNES, finite differences, Hessian, coloring, sparse

173: .seealso: SNESSetHessian()
174: @*/
175: int SNESDefaultComputeHessianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
176: {

180:   SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,ctx);
181:   return(0);
182: }