xref: /petsc/src/snes/interface/snesj.c (revision 84a7bf8dfd47cdf3a306d6ad5f80ccc190e1f109)
111320018SBarry Smith 
211320018SBarry Smith #ifndef lint
3*84a7bf8dSLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.35 1996/09/25 02:47:13 curfman Exp curfman $";
411320018SBarry Smith #endif
511320018SBarry Smith 
6b1f0a012SBarry Smith #include "draw.h"    /*I  "draw.h"  I*/
770f55243SBarry Smith #include "src/snes/snesimpl.h"    /*I  "snes.h"  I*/
811320018SBarry Smith 
94b828684SBarry Smith /*@C
10*84a7bf8dSLois 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:
1723242f5aSBarry Smith .  J - Jacobian
1823242f5aSBarry Smith .  B - preconditioner, same as Jacobian
19da88328cSLois Curfman McInnes .  flag - matrix flag
2011320018SBarry Smith 
21ad960d00SLois Curfman McInnes    Options Database Key:
22ad960d00SLois Curfman McInnes $  -snes_fd
23ad960d00SLois Curfman McInnes 
245f3c43d9SLois Curfman McInnes    Notes:
255f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
265f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
275f3c43d9SLois Curfman McInnes    SNESDefaultComputeJacobian() is not recommended for general use
285f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
295f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3011320018SBarry Smith 
315f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian
325f3c43d9SLois Curfman McInnes 
335f3c43d9SLois Curfman McInnes .seealso: SNESSetJacobian(), SNESTestJacobian()
3411320018SBarry Smith @*/
355334005bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
3611320018SBarry Smith {
3723242f5aSBarry Smith   Vec      j1,j2,x2;
3823242f5aSBarry Smith   int      i,ierr,N,start,end,j;
39bbb6d6a8SBarry Smith   Scalar   dx, mone = -1.0,*y,scale,*xx,wscale;
405334005bSBarry Smith   double   amax, epsilon = 1.e-8; /* assumes double precision */
4154d73c9fSLois Curfman McInnes   double   dx_min = 1.e-16, dx_par = 1.e-1;
42bbb6d6a8SBarry Smith   MPI_Comm comm;
430521c3abSLois Curfman McInnes   int      (*eval_fct)(SNES,Vec,Vec);
440521c3abSLois Curfman McInnes 
45a86d99e1SLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS)
460521c3abSLois Curfman McInnes     eval_fct = SNESComputeFunction;
47a86d99e1SLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
480521c3abSLois Curfman McInnes     eval_fct = SNESComputeGradient;
490521c3abSLois Curfman McInnes   else SETERRQ(1,"SNESDefaultComputeJacobian: Invalid method class");
5023242f5aSBarry Smith 
51bbb6d6a8SBarry Smith   PetscObjectGetComm((PetscObject)x1,&comm);
52336a5e98SBarry Smith   MatZeroEntries(*J);
53aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
54aa79bc6dSLois Curfman McInnes     ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr);
55aa79bc6dSLois Curfman McInnes     snes->nvwork = 3;
56aa79bc6dSLois Curfman McInnes     PLogObjectParents(snes,3,snes->vwork);
57aa79bc6dSLois Curfman McInnes   }
58aa79bc6dSLois Curfman McInnes   j1 = snes->vwork[0]; j2 = snes->vwork[1]; x2 = snes->vwork[2];
5923242f5aSBarry Smith 
6078b31e54SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
6178b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
6239e2f89bSBarry Smith   VecGetArray(x1,&xx);
630521c3abSLois Curfman McInnes   ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr);
64c005e166SLois Curfman McInnes 
65c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
66c005e166SLois Curfman McInnes       x1 = current iterate, j1 = F(x1)
67c005e166SLois Curfman McInnes       x2 = perturbed iterate, j2 = F(x2)
68c005e166SLois Curfman McInnes    */
6939e2f89bSBarry Smith   for ( i=0; i<N; i++ ) {
7078b31e54SBarry Smith     ierr = VecCopy(x1,x2); CHKERRQ(ierr);
7123242f5aSBarry Smith     if ( i>= start && i<end) {
7239e2f89bSBarry Smith       dx = xx[i-start];
7319a167f6SBarry Smith #if !defined(PETSC_COMPLEX)
7454d73c9fSLois Curfman McInnes       if (dx < dx_min && dx >= 0.0) dx = dx_par;
7554d73c9fSLois Curfman McInnes       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
7619a167f6SBarry Smith #else
7754d73c9fSLois Curfman McInnes       if (abs(dx) < dx_min && real(dx) >= 0.0) dx = dx_par;
7854d73c9fSLois Curfman McInnes       else if (real(dx) < 0.0 && abs(dx) < dx_min) dx = -dx_par;
7919a167f6SBarry Smith #endif
8039e2f89bSBarry Smith       dx *= epsilon;
8174f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
82dbb450caSBarry Smith       VecSetValues(x2,1,&i,&dx,ADD_VALUES);
8311320018SBarry Smith     }
84bbb6d6a8SBarry Smith     else {
85bbb6d6a8SBarry Smith       wscale = 0.0;
86bbb6d6a8SBarry Smith     }
870521c3abSLois Curfman McInnes     ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr);
8878b31e54SBarry Smith     ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr);
89c005e166SLois Curfman McInnes     /* Communicate scale to all processors */
90bbb6d6a8SBarry Smith #if !defined(PETSC_COMPLEX)
91bbb6d6a8SBarry Smith     MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);
92bbb6d6a8SBarry Smith #else
93bbb6d6a8SBarry Smith     MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);
94bbb6d6a8SBarry Smith #endif
9523242f5aSBarry Smith     VecScale(&scale,j2);
9623242f5aSBarry Smith     VecGetArray(j2,&y);
97cddf8d76SBarry Smith     VecNorm(j2,NORM_INFINITY,&amax); amax *= 1.e-14;
9823242f5aSBarry Smith     for ( j=start; j<end; j++ ) {
99cddf8d76SBarry Smith       if (PetscAbsScalar(y[j-start]) > amax) {
100dbb450caSBarry Smith         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr);
10123242f5aSBarry Smith       }
10223242f5aSBarry Smith     }
10323242f5aSBarry Smith     VecRestoreArray(j2,&y);
10423242f5aSBarry Smith   }
1050989deb7SLois Curfman McInnes   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1060989deb7SLois Curfman McInnes   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
107595b82d2SBarry Smith   *flag =  DIFFERENT_NONZERO_PATTERN;
10811320018SBarry Smith   return 0;
10911320018SBarry Smith }
11011320018SBarry Smith 
1110521c3abSLois Curfman McInnes /*@C
112*84a7bf8dSLois Curfman McInnes    SNESDefaultComputeHessian - Computes the Hessian using finite differences.
1130521c3abSLois Curfman McInnes 
1140521c3abSLois Curfman McInnes    Input Parameters:
1150521c3abSLois Curfman McInnes .  x1 - compute Hessian at this point
1160521c3abSLois Curfman McInnes .  ctx - application's gradient context, as set with SNESSetGradient()
1170521c3abSLois Curfman McInnes 
1180521c3abSLois Curfman McInnes    Output Parameters:
1190521c3abSLois Curfman McInnes .  J - Hessian
1200521c3abSLois Curfman McInnes .  B - preconditioner, same as Hessian
1210521c3abSLois Curfman McInnes .  flag - matrix flag
1220521c3abSLois Curfman McInnes 
1230521c3abSLois Curfman McInnes    Options Database Key:
1240521c3abSLois Curfman McInnes $  -snes_fd
1250521c3abSLois Curfman McInnes 
1260521c3abSLois Curfman McInnes    Notes:
1270521c3abSLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
1280521c3abSLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
1290521c3abSLois Curfman McInnes    SNESDefaultComputeHessian() is not recommended for general use
1300521c3abSLois Curfman McInnes    in large-scale applications, It can be useful in checking the
1310521c3abSLois Curfman McInnes    correctness of a user-provided Hessian.
1320521c3abSLois Curfman McInnes 
1330521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian
1340521c3abSLois Curfman McInnes 
1350521c3abSLois Curfman McInnes .seealso: SNESSetHessian(), SNESTestHessian()
1360521c3abSLois Curfman McInnes @*/
1370521c3abSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
1380521c3abSLois Curfman McInnes {
1390521c3abSLois Curfman McInnes   return SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);
1400521c3abSLois Curfman McInnes }
141