xref: /petsc/src/snes/interface/snesj.c (revision 336a5e98162d189faf560c36b5de145b11386e01)
111320018SBarry Smith 
239e2f89bSBarry Smith 
311320018SBarry Smith #ifndef lint
4*336a5e98SBarry Smith static char vcid[] = "$Id: snesj.c,v 1.7 1995/05/05 03:51:29 bsmith Exp bsmith $";
511320018SBarry Smith #endif
611320018SBarry Smith 
711320018SBarry Smith #include "draw.h"
839e2f89bSBarry Smith #include "snes.h"
911320018SBarry Smith #include "options.h"
1011320018SBarry Smith 
1111320018SBarry Smith /*@
1223242f5aSBarry Smith    SNESDefaultComputeJacobian - Computes Jacobian using finite
1323242f5aSBarry Smith        differences. Slow and expensive.
1411320018SBarry Smith 
1511320018SBarry Smith  Input Parameters:
1623242f5aSBarry Smith .  x - compute Jacobian at this point
1723242f5aSBarry Smith .  ctx - applications Function context
1811320018SBarry Smith 
1911320018SBarry Smith   Output Parameters:
2023242f5aSBarry Smith .  J - Jacobian
2123242f5aSBarry Smith .  B - preconditioner, same as Jacobian
2211320018SBarry Smith 
2323242f5aSBarry Smith .keywords: finite differences, Jacobian
2411320018SBarry Smith 
2523242f5aSBarry Smith .seealso: SNESSetJacobian, SNESTestJacobian
2611320018SBarry Smith @*/
2723242f5aSBarry Smith int SNESDefaultComputeJacobian(SNES snes, Vec x1,Mat *J,Mat *B,int *flag,
2823242f5aSBarry Smith                                void *ctx)
2911320018SBarry Smith {
3023242f5aSBarry Smith   Vec    j1,j2,x2;
3123242f5aSBarry Smith   int    i,ierr,N,start,end,j;
3239e2f89bSBarry Smith   Scalar dx, mone = -1.0,*y,scale,*xx;
33*336a5e98SBarry Smith   double epsilon = 1.e-8,amax; /* assumes double precision */
3423242f5aSBarry Smith 
35*336a5e98SBarry Smith   MatZeroEntries(*J);
366469c4f9SBarry Smith   ierr = VecDuplicate(x1,&j1); CHKERR(ierr);
376469c4f9SBarry Smith   ierr = VecDuplicate(x1,&j2); CHKERR(ierr);
386469c4f9SBarry Smith   ierr = VecDuplicate(x1,&x2); CHKERR(ierr);
3923242f5aSBarry Smith 
4023242f5aSBarry Smith   ierr = VecGetSize(x1,&N); CHKERR(ierr);
4123242f5aSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERR(ierr);
4239e2f89bSBarry Smith   VecGetArray(x1,&xx);
4323242f5aSBarry Smith   ierr = SNESComputeFunction(snes,x1,j1); CHKERR(ierr);
4439e2f89bSBarry Smith   for ( i=0; i<N; i++ ) {
45*336a5e98SBarry Smith     ierr = VecCopy(x1,x2); CHKERR(ierr);
4623242f5aSBarry Smith     if ( i>= start && i<end) {
4739e2f89bSBarry Smith       dx = xx[i-start];
48*336a5e98SBarry Smith       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
49*336a5e98SBarry Smith       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
5039e2f89bSBarry Smith       dx *= epsilon;
5139e2f89bSBarry Smith       scale = -1.0/dx;
529d00d63dSBarry Smith       VecSetValues(x2,1,&i,&dx,ADDVALUES);
5311320018SBarry Smith     }
5423242f5aSBarry Smith     ierr = SNESComputeFunction(snes,x2,j2); CHKERR(ierr);
5523242f5aSBarry Smith     ierr = VecAXPY(&mone,j1,j2); CHKERR(ierr);
5623242f5aSBarry Smith     VecScale(&scale,j2);
5723242f5aSBarry Smith     VecGetArray(j2,&y);
58*336a5e98SBarry Smith     VecAMax(j2,0,&amax); amax *= 1.e-14;
5923242f5aSBarry Smith     for ( j=start; j<end; j++ ) {
60*336a5e98SBarry Smith       if (y[j-start] > amax || y[j-start] < -amax) {
61edae2e7dSBarry Smith         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERR(ierr);
6223242f5aSBarry Smith       }
6323242f5aSBarry Smith     }
6423242f5aSBarry Smith     VecRestoreArray(j2,&y);
6523242f5aSBarry Smith   }
66ee50ffe9SBarry Smith   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
6723242f5aSBarry Smith   VecDestroy(x2); VecDestroy(j1); VecDestroy(j2);
68ee50ffe9SBarry Smith   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
6911320018SBarry Smith   return 0;
7011320018SBarry Smith }
7111320018SBarry Smith 
72