Actual source code: snestest.c
1: #define PETSCSNES_DLL
3: #include private/snesimpl.h
5: typedef struct {
6: PetscTruth complete_print;
7: } SNES_Test;
9: /*
10: SNESSolve_Test - Tests whether a hand computed Jacobian
11: matches one compute via finite differences.
12: */
15: PetscErrorCode SNESSolve_Test(SNES snes)
16: {
17: Mat A = snes->jacobian,B;
18: Vec x = snes->vec_sol, f = snes->vec_func;
20: PetscInt i;
21: MatStructure flg;
22: PetscReal nrm,gnorm;
23: SNES_Test *neP = (SNES_Test*)snes->data;
27: if (A != snes->jacobian_pre) {
28: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot test with alternative preconditioner");
29: }
31: PetscPrintf(((PetscObject)snes)->comm,"Testing hand-coded Jacobian, if the ratio is\n");
32: PetscPrintf(((PetscObject)snes)->comm,"O(1.e-8), the hand-coded Jacobian is probably correct.\n");
33: if (!neP->complete_print) {
34: PetscPrintf(((PetscObject)snes)->comm,"Run with -snes_test_display to show difference\n");
35: PetscPrintf(((PetscObject)snes)->comm,"of hand-coded and finite difference Jacobian.\n");
36: }
38: for (i=0; i<3; i++) {
39: if (i == 1) {VecSet(x,-1.0);}
40: else if (i == 2) {VecSet(x,1.0);}
41:
42: /* evaluate the function at this point because SNESDefaultComputeJacobianColor() assumes that the function has been evaluated and put into snes->vec_func */
43: SNESComputeFunction(snes,x,f);
45: /* compute both versions of Jacobian */
46: SNESComputeJacobian(snes,x,&A,&A,&flg);
47: if (!i) {MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&B);}
48: SNESDefaultComputeJacobian(snes,x,&B,&B,&flg,snes->funP);
49: if (neP->complete_print) {
50: MPI_Comm comm;
51: PetscViewer viewer;
52: PetscPrintf(((PetscObject)snes)->comm,"Finite difference Jacobian\n");
53: PetscObjectGetComm((PetscObject)B,&comm);
54: PetscViewerASCIIGetStdout(comm,&viewer);
55: MatView(B,viewer);
56: }
57: /* compare */
58: MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
59: MatNorm(B,NORM_FROBENIUS,&nrm);
60: MatNorm(A,NORM_FROBENIUS,&gnorm);
61: if (neP->complete_print) {
62: MPI_Comm comm;
63: PetscViewer viewer;
64: PetscPrintf(((PetscObject)snes)->comm,"Hand-coded Jacobian\n");
65: PetscObjectGetComm((PetscObject)B,&comm);
66: PetscViewerASCIIGetStdout(comm,&viewer);
67: MatView(A,viewer);
68: }
69: if (!gnorm) gnorm = 1; /* just in case */
70: PetscPrintf(((PetscObject)snes)->comm,"Norm of matrix ratio %G difference %G\n",nrm/gnorm,nrm);
71: }
72: MatDestroy(B);
73: /*
74: Return error code cause Jacobian not good
75: */
76: PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
77: }
78: /* ------------------------------------------------------------ */
81: PetscErrorCode SNESDestroy_Test(SNES snes)
82: {
85: PetscFree(snes->data);
86: return(0);
87: }
91: static PetscErrorCode SNESSetFromOptions_Test(SNES snes)
92: {
93: SNES_Test *ls = (SNES_Test *)snes->data;
98: PetscOptionsHead("Hand-coded Jacobian tester options");
99: PetscOptionsName("-snes_test_display","Display difference between approximate and handcoded Jacobian","None",&ls->complete_print);
100: PetscOptionsTail();
101: return(0);
102: }
104: /* ------------------------------------------------------------ */
105: /*MC
106: SNESTEST - Test hand-coded Jacobian against finite difference Jacobian
108: Options Database:
109: . -snes_test_display Display difference between approximate and hand-coded Jacobian
111: Level: intermediate
113: .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
115: M*/
119: PetscErrorCode SNESCreate_Test(SNES snes)
120: {
121: SNES_Test *neP;
125: snes->ops->setup = 0;
126: snes->ops->solve = SNESSolve_Test;
127: snes->ops->destroy = SNESDestroy_Test;
128: snes->ops->setfromoptions = SNESSetFromOptions_Test;
129: snes->ops->view = 0;
131: ierr = PetscNewLog(snes,SNES_Test,&neP);
132: snes->data = (void*)neP;
133: neP->complete_print = PETSC_FALSE;
134: return(0);
135: }