xref: /petsc/src/snes/interface/snesj.c (revision 88c956ad8fd32d7b74ab816df739cc596c57dec4)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*88c956adSLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.48 1998/04/13 17:55:33 bsmith Exp curfman $";
311320018SBarry Smith #endif
411320018SBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"    /*I  "snes.h"  I*/
611320018SBarry Smith 
75615d1e5SSatish Balay #undef __FUNC__
85615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeJacobian"
94b828684SBarry Smith /*@C
1084a7bf8dSLois Curfman McInnes    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
1111320018SBarry Smith 
1211320018SBarry Smith    Input Parameters:
13da88328cSLois Curfman McInnes .  x1 - compute Jacobian at this point
14da88328cSLois Curfman McInnes .  ctx - application's function context, as set with SNESSetFunction()
1511320018SBarry Smith 
1611320018SBarry Smith    Output Parameters:
17b4fc646aSLois Curfman McInnes .  J - Jacobian matrix (not altered in this routine)
18b4fc646aSLois Curfman McInnes .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
19b4fc646aSLois Curfman McInnes .  flag - flag indicating whether the matrix sparsity structure has changed
2011320018SBarry Smith 
21fee21e36SBarry Smith    Collective on SNES
22fee21e36SBarry Smith 
23ad960d00SLois Curfman McInnes    Options Database Key:
24ad960d00SLois Curfman McInnes $  -snes_fd
25ad960d00SLois Curfman McInnes 
265f3c43d9SLois Curfman McInnes    Notes:
275f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
285f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
295f3c43d9SLois Curfman McInnes    SNESDefaultComputeJacobian() is not recommended for general use
305f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
315f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3211320018SBarry Smith 
33b4fc646aSLois Curfman McInnes    An alternative routine that uses coloring to explot matrix sparsity is
34b4fc646aSLois Curfman McInnes    SNESDefaultComputeJacobianWithColoring().
35b4fc646aSLois Curfman McInnes 
365f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian
375f3c43d9SLois Curfman McInnes 
38b4fc646aSLois Curfman McInnes .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianWithColoring()
3911320018SBarry Smith @*/
405334005bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
4111320018SBarry Smith {
42*88c956adSLois Curfman McInnes   Vec      j1a,j2a,x2;
4323242f5aSBarry Smith   int      i,ierr,N,start,end,j;
44bbb6d6a8SBarry Smith   Scalar   dx, mone = -1.0,*y,scale,*xx,wscale;
455334005bSBarry Smith   double   amax, epsilon = 1.e-8; /* assumes double precision */
4654d73c9fSLois Curfman McInnes   double   dx_min = 1.e-16, dx_par = 1.e-1;
47bbb6d6a8SBarry Smith   MPI_Comm comm;
480521c3abSLois Curfman McInnes   int      (*eval_fct)(SNES,Vec,Vec);
490521c3abSLois Curfman McInnes 
503a40ed3dSBarry Smith   PetscFunctionBegin;
513a40ed3dSBarry Smith   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction;
523a40ed3dSBarry Smith   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient;
53a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
5423242f5aSBarry Smith 
55bbb6d6a8SBarry Smith   PetscObjectGetComm((PetscObject)x1,&comm);
56b4fc646aSLois Curfman McInnes   MatZeroEntries(*B);
57aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
58aa79bc6dSLois Curfman McInnes     ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr);
59aa79bc6dSLois Curfman McInnes     snes->nvwork = 3;
60aa79bc6dSLois Curfman McInnes     PLogObjectParents(snes,3,snes->vwork);
61aa79bc6dSLois Curfman McInnes   }
62*88c956adSLois Curfman McInnes   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
6323242f5aSBarry Smith 
6478b31e54SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
6578b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
6639e2f89bSBarry Smith   VecGetArray(x1,&xx);
67*88c956adSLois Curfman McInnes   ierr = eval_fct(snes,x1,j1a); CHKERRQ(ierr);
68c005e166SLois Curfman McInnes 
69c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
70*88c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
71*88c956adSLois 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) {
7639e2f89bSBarry Smith       dx = xx[i-start];
773a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
7854d73c9fSLois Curfman McInnes       if (dx < dx_min && dx >= 0.0) dx = dx_par;
7954d73c9fSLois Curfman McInnes       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
8019a167f6SBarry Smith #else
8154d73c9fSLois Curfman McInnes       if (abs(dx) < dx_min && real(dx) >= 0.0) dx = dx_par;
8254d73c9fSLois Curfman McInnes       else if (real(dx) < 0.0 && abs(dx) < dx_min) dx = -dx_par;
8319a167f6SBarry Smith #endif
8439e2f89bSBarry Smith       dx *= epsilon;
8574f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
86dbb450caSBarry Smith       VecSetValues(x2,1,&i,&dx,ADD_VALUES);
8711320018SBarry Smith     }
88bbb6d6a8SBarry Smith     else {
89bbb6d6a8SBarry Smith       wscale = 0.0;
90bbb6d6a8SBarry Smith     }
91*88c956adSLois Curfman McInnes     ierr = eval_fct(snes,x2,j2a); CHKERRQ(ierr);
92*88c956adSLois Curfman McInnes     ierr = VecAXPY(&mone,j1a,j2a); CHKERRQ(ierr);
93c005e166SLois Curfman McInnes     /* Communicate scale to all processors */
943a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
95ca161407SBarry Smith     ierr = MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
96bbb6d6a8SBarry Smith #else
97ca161407SBarry Smith     ierr = MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
98bbb6d6a8SBarry Smith #endif
99*88c956adSLois Curfman McInnes     VecScale(&scale,j2a);
100*88c956adSLois Curfman McInnes     VecGetArray(j2a,&y);
101*88c956adSLois Curfman McInnes     VecNorm(j2a,NORM_INFINITY,&amax); amax *= 1.e-14;
10223242f5aSBarry Smith     for ( j=start; j<end; j++ ) {
103cddf8d76SBarry Smith       if (PetscAbsScalar(y[j-start]) > amax) {
104b4fc646aSLois Curfman McInnes         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr);
10523242f5aSBarry Smith       }
10623242f5aSBarry Smith     }
107*88c956adSLois Curfman McInnes     VecRestoreArray(j2a,&y);
10823242f5aSBarry Smith   }
109b4fc646aSLois Curfman McInnes   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
110b4fc646aSLois Curfman McInnes   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
111595b82d2SBarry Smith   *flag =  DIFFERENT_NONZERO_PATTERN;
1123a40ed3dSBarry Smith   PetscFunctionReturn(0);
11311320018SBarry Smith }
11411320018SBarry Smith 
1155615d1e5SSatish Balay #undef __FUNC__
1165615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeHessian"
1170521c3abSLois Curfman McInnes /*@C
11884a7bf8dSLois Curfman McInnes    SNESDefaultComputeHessian - Computes the Hessian using finite differences.
1190521c3abSLois Curfman McInnes 
1200521c3abSLois Curfman McInnes    Input Parameters:
1210521c3abSLois Curfman McInnes .  x1 - compute Hessian at this point
1220521c3abSLois Curfman McInnes .  ctx - application's gradient context, as set with SNESSetGradient()
1230521c3abSLois Curfman McInnes 
1240521c3abSLois Curfman McInnes    Output Parameters:
125b4fc646aSLois Curfman McInnes .  J - Hessian matrix (not altered in this routine)
126b4fc646aSLois Curfman McInnes .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
127b4fc646aSLois Curfman McInnes .  flag - flag indicating whether the matrix sparsity structure has changed
1280521c3abSLois Curfman McInnes 
129fee21e36SBarry Smith    Collective on SNES
130fee21e36SBarry Smith 
1310521c3abSLois Curfman McInnes    Options Database Key:
1320521c3abSLois Curfman McInnes $  -snes_fd
1330521c3abSLois Curfman McInnes 
1340521c3abSLois Curfman McInnes    Notes:
1350521c3abSLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
1360521c3abSLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
1370521c3abSLois Curfman McInnes    SNESDefaultComputeHessian() is not recommended for general use
1380521c3abSLois Curfman McInnes    in large-scale applications, It can be useful in checking the
1390521c3abSLois Curfman McInnes    correctness of a user-provided Hessian.
1400521c3abSLois Curfman McInnes 
1410521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian
1420521c3abSLois Curfman McInnes 
1439984d6fbSBarry Smith .seealso: SNESSetHessian()
1440521c3abSLois Curfman McInnes @*/
1450521c3abSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
1460521c3abSLois Curfman McInnes {
1473a40ed3dSBarry Smith   int ierr;
1483a40ed3dSBarry Smith 
1493a40ed3dSBarry Smith   PetscFunctionBegin;
1503a40ed3dSBarry Smith   ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
1513a40ed3dSBarry Smith   PetscFunctionReturn(0);
1520521c3abSLois Curfman McInnes }
153