xref: /petsc/src/snes/interface/snesj.c (revision ec12a082eca8625a8f0fdf4e67e9cce9e40409bf)
111320018SBarry Smith 
211320018SBarry Smith #ifndef lint
3*ec12a082SBarry Smith static char vcid[] = "$Id: snesj.c,v 1.20 1995/08/14 23:12:10 curfman Exp bsmith $";
411320018SBarry Smith #endif
511320018SBarry Smith 
6b1f0a012SBarry Smith #include "draw.h"    /*I  "draw.h"  I*/
7b1f0a012SBarry Smith #include "snes.h"    /*I  "snes.h"  I*/
811320018SBarry Smith 
911320018SBarry Smith /*@
105f3c43d9SLois Curfman McInnes    SNESDefaultComputeJacobian - Computes the Jacobian using finite
115f3c43d9SLois Curfman McInnes    differences.
1211320018SBarry Smith 
1311320018SBarry Smith    Input Parameters:
14da88328cSLois Curfman McInnes .  x1 - compute Jacobian at this point
15da88328cSLois Curfman McInnes .  ctx - application's function context, as set with SNESSetFunction()
1611320018SBarry Smith 
1711320018SBarry Smith    Output Parameters:
1823242f5aSBarry Smith .  J - Jacobian
1923242f5aSBarry Smith .  B - preconditioner, same as Jacobian
20da88328cSLois Curfman McInnes .  flag - matrix flag
2111320018SBarry Smith 
22ad960d00SLois Curfman McInnes    Options Database Key:
23ad960d00SLois Curfman McInnes $  -snes_fd
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 
325f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian
335f3c43d9SLois Curfman McInnes 
345f3c43d9SLois Curfman McInnes .seealso: SNESSetJacobian(), SNESTestJacobian()
3511320018SBarry Smith @*/
36eccfb7ebSLois Curfman McInnes int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,
37eccfb7ebSLois Curfman McInnes                                MatStructure *flag,void *ctx)
3811320018SBarry Smith {
3923242f5aSBarry Smith   Vec      j1,j2,x2;
4023242f5aSBarry Smith   int      i,ierr,N,start,end,j;
41bbb6d6a8SBarry Smith   Scalar   dx, mone = -1.0,*y,scale,*xx,wscale;
42336a5e98SBarry Smith   double   epsilon = 1.e-8,amax; /* assumes double precision */
43bbb6d6a8SBarry Smith   MPI_Comm comm;
4423242f5aSBarry Smith 
45bbb6d6a8SBarry Smith   PetscObjectGetComm((PetscObject)x1,&comm);
46336a5e98SBarry Smith   MatZeroEntries(*J);
4778b31e54SBarry Smith   ierr = VecDuplicate(x1,&j1); CHKERRQ(ierr);
4878b31e54SBarry Smith   ierr = VecDuplicate(x1,&j2); CHKERRQ(ierr);
4978b31e54SBarry Smith   ierr = VecDuplicate(x1,&x2); CHKERRQ(ierr);
50d370d78aSBarry Smith   PLogObjectParent(snes,j1); PLogObjectParent(snes,j2);
51393d2d9aSLois Curfman McInnes   PLogObjectParent(snes,x2);
5223242f5aSBarry Smith 
5378b31e54SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
5478b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
5539e2f89bSBarry Smith   VecGetArray(x1,&xx);
5678b31e54SBarry Smith   ierr = SNESComputeFunction(snes,x1,j1); CHKERRQ(ierr);
5739e2f89bSBarry Smith   for ( i=0; i<N; i++ ) {
5878b31e54SBarry Smith     ierr = VecCopy(x1,x2); CHKERRQ(ierr);
5923242f5aSBarry Smith     if ( i>= start && i<end) {
6039e2f89bSBarry Smith       dx = xx[i-start];
6119a167f6SBarry Smith #if !defined(PETSC_COMPLEX)
62336a5e98SBarry Smith       if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1;
63336a5e98SBarry Smith       else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1;
6419a167f6SBarry Smith #else
6519a167f6SBarry Smith       if (abs(dx) < 1.e-16 && real(dx) >= 0.0) dx = 1.e-1;
66*ec12a082SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < 1.e-16) dx = -1.e-1;
6719a167f6SBarry Smith #endif
6839e2f89bSBarry Smith       dx *= epsilon;
69bbb6d6a8SBarry Smith       wscale = -1.0/dx;
709d00d63dSBarry Smith       VecSetValues(x2,1,&i,&dx,ADDVALUES);
7111320018SBarry Smith     }
72bbb6d6a8SBarry Smith     else {
73bbb6d6a8SBarry Smith       wscale = 0.0;
74bbb6d6a8SBarry Smith     }
7578b31e54SBarry Smith     ierr = SNESComputeFunction(snes,x2,j2); CHKERRQ(ierr);
7678b31e54SBarry Smith     ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr);
77bbb6d6a8SBarry Smith /* communicate scale to all processors */
78bbb6d6a8SBarry Smith #if !defined(PETSC_COMPLEX)
79bbb6d6a8SBarry Smith     MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);
80bbb6d6a8SBarry Smith #else
81bbb6d6a8SBarry Smith     MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);
82bbb6d6a8SBarry Smith #endif
8323242f5aSBarry Smith     VecScale(&scale,j2);
8423242f5aSBarry Smith     VecGetArray(j2,&y);
85336a5e98SBarry Smith     VecAMax(j2,0,&amax); amax *= 1.e-14;
8623242f5aSBarry Smith     for ( j=start; j<end; j++ ) {
8719a167f6SBarry Smith #if defined(PETSC_COMPLEX)
8819a167f6SBarry Smith       if (abs(y[j-start]) > amax) {
8919a167f6SBarry Smith #else
90336a5e98SBarry Smith       if (y[j-start] > amax || y[j-start] < -amax) {
9119a167f6SBarry Smith #endif
9278b31e54SBarry Smith         ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERRQ(ierr);
9323242f5aSBarry Smith       }
9423242f5aSBarry Smith     }
9523242f5aSBarry Smith     VecRestoreArray(j2,&y);
9623242f5aSBarry Smith   }
97ee50ffe9SBarry Smith   MatAssemblyBegin(*J,FINAL_ASSEMBLY);
9823242f5aSBarry Smith   VecDestroy(x2); VecDestroy(j1); VecDestroy(j2);
99ee50ffe9SBarry Smith   MatAssemblyEnd(*J,FINAL_ASSEMBLY);
10011320018SBarry Smith   return 0;
10111320018SBarry Smith }
10211320018SBarry Smith 
103