1 2 3 #ifndef lint 4 static char vcid[] = "$Id: snesj.c,v 1.9 1995/05/11 17:53:49 curfman Exp bsmith $"; 5 #endif 6 7 #include "draw.h" 8 #include "snes.h" 9 #include "options.h" 10 11 /*@ 12 SNESDefaultComputeJacobian - Computes the Jacobian using finite 13 differences. 14 15 Input Parameters: 16 . x - compute Jacobian at this point 17 . ctx - application's Function context 18 19 Output Parameters: 20 . J - Jacobian 21 . B - preconditioner, same as Jacobian 22 23 Notes: 24 This routine is slow and expensive, and is not currently optimized 25 to take advantage of sparsity in the problem. Although 26 SNESDefaultComputeJacobian() is not recommended for general use 27 in large-scale applications, It can be useful in checking the 28 correctness of a user-provided Jacobian. 29 30 .keywords: SNES, finite differences, Jacobian 31 32 .seealso: SNESSetJacobian(), SNESTestJacobian() 33 @*/ 34 int SNESDefaultComputeJacobian(SNES snes, Vec x1,Mat *J,Mat *B,int *flag, 35 void *ctx) 36 { 37 Vec j1,j2,x2; 38 int i,ierr,N,start,end,j; 39 Scalar dx, mone = -1.0,*y,scale,*xx; 40 double epsilon = 1.e-8,amax; /* assumes double precision */ 41 42 MatZeroEntries(*J); 43 ierr = VecDuplicate(x1,&j1); CHKERR(ierr); 44 ierr = VecDuplicate(x1,&j2); CHKERR(ierr); 45 ierr = VecDuplicate(x1,&x2); CHKERR(ierr); 46 47 ierr = VecGetSize(x1,&N); CHKERR(ierr); 48 ierr = VecGetOwnershipRange(x1,&start,&end); CHKERR(ierr); 49 VecGetArray(x1,&xx); 50 ierr = SNESComputeFunction(snes,x1,j1); CHKERR(ierr); 51 for ( i=0; i<N; i++ ) { 52 ierr = VecCopy(x1,x2); CHKERR(ierr); 53 if ( i>= start && i<end) { 54 dx = xx[i-start]; 55 if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1; 56 else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1; 57 dx *= epsilon; 58 scale = -1.0/dx; 59 VecSetValues(x2,1,&i,&dx,ADDVALUES); 60 } 61 ierr = SNESComputeFunction(snes,x2,j2); CHKERR(ierr); 62 ierr = VecAXPY(&mone,j1,j2); CHKERR(ierr); 63 VecScale(&scale,j2); 64 VecGetArray(j2,&y); 65 VecAMax(j2,0,&amax); amax *= 1.e-14; 66 for ( j=start; j<end; j++ ) { 67 if (y[j-start] > amax || y[j-start] < -amax) { 68 ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERTVALUES); CHKERR(ierr); 69 } 70 } 71 VecRestoreArray(j2,&y); 72 } 73 MatAssemblyBegin(*J,FINAL_ASSEMBLY); 74 VecDestroy(x2); VecDestroy(j1); VecDestroy(j2); 75 MatAssemblyEnd(*J,FINAL_ASSEMBLY); 76 return 0; 77 } 78 79