xref: /petsc/src/snes/interface/snesj.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1*29bbc08cSBarry Smith /*$Id: snesj.c,v 1.67 2000/08/31 17:10:59 bsmith Exp bsmith $*/
211320018SBarry Smith 
3e090d566SSatish Balay #include "src/snes/snesimpl.h"    /*I  "petscsnes.h"  I*/
411320018SBarry Smith 
55615d1e5SSatish Balay #undef __FUNC__
6b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDefaultComputeJacobian"
74b828684SBarry Smith /*@C
884a7bf8dSLois Curfman McInnes    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
911320018SBarry Smith 
10fee21e36SBarry Smith    Collective on SNES
11fee21e36SBarry Smith 
12c7afd0dbSLois Curfman McInnes    Input Parameters:
13c7afd0dbSLois Curfman McInnes +  x1 - compute Jacobian at this point
14c7afd0dbSLois Curfman McInnes -  ctx - application's function context, as set with SNESSetFunction()
15c7afd0dbSLois Curfman McInnes 
16c7afd0dbSLois Curfman McInnes    Output Parameters:
17c7afd0dbSLois Curfman McInnes +  J - Jacobian matrix (not altered in this routine)
18c7afd0dbSLois Curfman McInnes .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
19c7afd0dbSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
20c7afd0dbSLois Curfman McInnes 
21ad960d00SLois Curfman McInnes    Options Database Key:
2208f75eefSBarry Smith +  -snes_fd - Activates SNESDefaultComputeJacobian()
2308f75eefSBarry Smith -  -snes_test_err - Square root of function error tolerance, default 1.e-8
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 @*/
412d0c0e3bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
4211320018SBarry Smith {
4388c956adSLois Curfman McInnes   Vec       j1a,j2a,x2;
4423242f5aSBarry Smith   int       i,ierr,N,start,end,j;
45bbb6d6a8SBarry Smith   Scalar    dx,mone = -1.0,*y,scale,*xx,wscale;
46329f5518SBarry Smith   PetscReal amax,epsilon = 1.e-8; /* assumes PetscReal precision */
47329f5518SBarry Smith   PetscReal dx_min = 1.e-16,dx_par = 1.e-1;
48bbb6d6a8SBarry Smith   MPI_Comm  comm;
49a305c92eSSatish Balay   int      (*eval_fct)(SNES,Vec,Vec)=0;
500521c3abSLois Curfman McInnes 
513a40ed3dSBarry Smith   PetscFunctionBegin;
5208f75eefSBarry Smith   ierr = OptionsGetDouble(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr);
533a40ed3dSBarry Smith   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction;
543a40ed3dSBarry Smith   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient;
55*29bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid method class");
5623242f5aSBarry Smith 
572d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
582d0c0e3bSBarry Smith   ierr = MatZeroEntries(*B);CHKERRQ(ierr);
59aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
60aa79bc6dSLois Curfman McInnes     ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr);
61aa79bc6dSLois Curfman McInnes     snes->nvwork = 3;
62aa79bc6dSLois Curfman McInnes     PLogObjectParents(snes,3,snes->vwork);
63aa79bc6dSLois Curfman McInnes   }
6488c956adSLois Curfman McInnes   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
6523242f5aSBarry Smith 
6678b31e54SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
6778b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
6888c956adSLois Curfman McInnes   ierr = eval_fct(snes,x1,j1a);CHKERRQ(ierr);
69c005e166SLois Curfman McInnes 
70c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
7188c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
7288c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
73c005e166SLois Curfman McInnes    */
7439e2f89bSBarry Smith   for (i=0; i<N; i++) {
7578b31e54SBarry Smith     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
7623242f5aSBarry Smith     if (i>= start && i<end) {
772e8a6d31SBarry Smith       ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
7839e2f89bSBarry Smith       dx = xx[i-start];
792e8a6d31SBarry Smith       ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
80aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8154d73c9fSLois Curfman McInnes       if (dx < dx_min && dx >= 0.0) dx = dx_par;
8254d73c9fSLois Curfman McInnes       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
8319a167f6SBarry Smith #else
84329f5518SBarry Smith       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
85329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
8619a167f6SBarry Smith #endif
8739e2f89bSBarry Smith       dx *= epsilon;
8874f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
892e8a6d31SBarry Smith       ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
906c783aadSBarry Smith     } else {
91bbb6d6a8SBarry Smith       wscale = 0.0;
92bbb6d6a8SBarry Smith     }
9388c956adSLois Curfman McInnes     ierr = eval_fct(snes,x2,j2a);CHKERRQ(ierr);
9488c956adSLois Curfman McInnes     ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr);
95c005e166SLois Curfman McInnes     /* Communicate scale to all processors */
966831982aSBarry Smith     ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
972e8a6d31SBarry Smith     ierr = VecScale(&scale,j2a);CHKERRQ(ierr);
982e8a6d31SBarry Smith     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
992e8a6d31SBarry Smith     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
10023242f5aSBarry Smith     for (j=start; j<end; j++) {
101cddf8d76SBarry Smith       if (PetscAbsScalar(y[j-start]) > amax) {
102b4fc646aSLois Curfman McInnes         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
10323242f5aSBarry Smith       }
10423242f5aSBarry Smith     }
1052e8a6d31SBarry Smith     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
10623242f5aSBarry Smith   }
107b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109595b82d2SBarry Smith   *flag =  DIFFERENT_NONZERO_PATTERN;
1103a40ed3dSBarry Smith   PetscFunctionReturn(0);
11111320018SBarry Smith }
11211320018SBarry Smith 
1135615d1e5SSatish Balay #undef __FUNC__
114b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDefaultComputeHessian"
1150521c3abSLois Curfman McInnes /*@C
11684a7bf8dSLois Curfman McInnes    SNESDefaultComputeHessian - Computes the Hessian using finite differences.
1170521c3abSLois Curfman McInnes 
118fee21e36SBarry Smith    Collective on SNES
119fee21e36SBarry Smith 
120c7afd0dbSLois Curfman McInnes    Input Parameters:
121c7afd0dbSLois Curfman McInnes +  x1 - compute Hessian at this point
122c7afd0dbSLois Curfman McInnes -  ctx - application's gradient context, as set with SNESSetGradient()
123c7afd0dbSLois Curfman McInnes 
124c7afd0dbSLois Curfman McInnes    Output Parameters:
125c7afd0dbSLois Curfman McInnes +  J - Hessian matrix (not altered in this routine)
126c7afd0dbSLois Curfman McInnes .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
127c7afd0dbSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
128c7afd0dbSLois Curfman McInnes 
1290521c3abSLois Curfman McInnes    Options Database Key:
130c7afd0dbSLois Curfman McInnes $  -snes_fd - Activates SNESDefaultComputeHessian()
1310521c3abSLois Curfman McInnes 
13215091d37SBarry Smith 
13315091d37SBarry Smith    Level: intermediate
13415091d37SBarry Smith 
1350521c3abSLois Curfman McInnes    Notes:
1360521c3abSLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
1370521c3abSLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
1380521c3abSLois Curfman McInnes    SNESDefaultComputeHessian() is not recommended for general use
1390521c3abSLois Curfman McInnes    in large-scale applications, It can be useful in checking the
1400521c3abSLois Curfman McInnes    correctness of a user-provided Hessian.
1410521c3abSLois Curfman McInnes 
1420521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian
1430521c3abSLois Curfman McInnes 
1449984d6fbSBarry Smith .seealso: SNESSetHessian()
1450521c3abSLois Curfman McInnes @*/
146329f5518SBarry Smith int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
1470521c3abSLois Curfman McInnes {
1483a40ed3dSBarry Smith   int ierr;
1493a40ed3dSBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
1513a40ed3dSBarry Smith   ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
1530521c3abSLois Curfman McInnes }
154bd45824fSLois Curfman McInnes 
155bd45824fSLois Curfman McInnes #undef __FUNC__
156b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDefaultComputeHessianColor"
157bd45824fSLois Curfman McInnes /*@C
158bd45824fSLois Curfman McInnes    SNESDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
159bd45824fSLois Curfman McInnes 
160bd45824fSLois Curfman McInnes    Collective on SNES
161bd45824fSLois Curfman McInnes 
162bd45824fSLois Curfman McInnes    Input Parameters:
163bd45824fSLois Curfman McInnes +  x1 - compute Hessian at this point
164bd45824fSLois Curfman McInnes -  ctx - application's gradient context, as set with SNESSetGradient()
165bd45824fSLois Curfman McInnes 
166bd45824fSLois Curfman McInnes    Output Parameters:
167bd45824fSLois Curfman McInnes +  J - Hessian matrix (not altered in this routine)
168bd45824fSLois Curfman McInnes .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
169bd45824fSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
170bd45824fSLois Curfman McInnes 
171bd45824fSLois Curfman McInnes     Options Database Keys:
172bd45824fSLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor()
173bd45824fSLois Curfman McInnes 
174bd45824fSLois Curfman McInnes    Level: intermediate
175bd45824fSLois Curfman McInnes 
176bd45824fSLois Curfman McInnes  .keywords: SNES, finite differences, Hessian, coloring, sparse
177bd45824fSLois Curfman McInnes 
178bd45824fSLois Curfman McInnes .seealso: SNESSetHessian()
179bd45824fSLois Curfman McInnes @*/
180329f5518SBarry Smith int SNESDefaultComputeHessianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
181bd45824fSLois Curfman McInnes {
182bd45824fSLois Curfman McInnes   int ierr;
183bd45824fSLois Curfman McInnes 
184bd45824fSLois Curfman McInnes   PetscFunctionBegin;
185bd45824fSLois Curfman McInnes   ierr = SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
186bd45824fSLois Curfman McInnes   PetscFunctionReturn(0);
187bd45824fSLois Curfman McInnes }
188bd45824fSLois Curfman McInnes 
189