Actual source code: snestest.c
1: #define PETSCSNES_DLL
3: #include 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;
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(snes->comm,"Testing hand-coded Jacobian, if the ratio is\n");
32: PetscPrintf(snes->comm,"O(1.e-8), the hand-coded Jacobian is probably correct.\n");
33: if (!neP->complete_print) {
34: PetscPrintf(snes->comm,"Run with -snes_test_display to show difference\n");
35: PetscPrintf(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: /* compute both versions of Jacobian */
43: SNESComputeJacobian(snes,x,&A,&A,&flg);
44: if (!i) {MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&B);}
45: SNESDefaultComputeJacobian(snes,x,&B,&B,&flg,snes->funP);
46: if (neP->complete_print) {
47: MPI_Comm comm;
48: PetscViewer viewer;
49: PetscPrintf(snes->comm,"Finite difference Jacobian\n");
50: PetscObjectGetComm((PetscObject)B,&comm);
51: PetscViewerASCIIGetStdout(comm,&viewer);
52: MatView(B,viewer);
53: }
54: /* compare */
55: MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
56: MatNorm(B,NORM_FROBENIUS,&nrm);
57: MatNorm(A,NORM_FROBENIUS,&gnorm);
58: if (neP->complete_print) {
59: MPI_Comm comm;
60: PetscViewer viewer;
61: PetscPrintf(snes->comm,"Hand-coded Jacobian\n");
62: PetscObjectGetComm((PetscObject)B,&comm);
63: PetscViewerASCIIGetStdout(comm,&viewer);
64: MatView(A,viewer);
65: }
66: PetscPrintf(snes->comm,"Norm of matrix ratio %G difference %G\n",nrm/gnorm,nrm);
67: }
68: MatDestroy(B);
69: /*
70: Return error code cause Jacobian not good
71: */
72: PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
73: }
74: /* ------------------------------------------------------------ */
77: PetscErrorCode SNESDestroy_Test(SNES snes)
78: {
81: PetscFree(snes->data);
82: return(0);
83: }
87: static PetscErrorCode SNESSetFromOptions_Test(SNES snes)
88: {
89: SNES_Test *ls = (SNES_Test *)snes->data;
94: PetscOptionsHead("Hand-coded Jacobian tester options");
95: PetscOptionsName("-snes_test_display","Display difference between approximate and handcoded Jacobian","None",&ls->complete_print);
96: PetscOptionsTail();
97: return(0);
98: }
100: /* ------------------------------------------------------------ */
101: /*MC
102: SNESTEST - Test hand-coded Jacobian against finite difference Jacobian
104: Options Database:
105: . -snes_test_display Display difference between approximate and hand-coded Jacobian
107: Level: intermediate
109: .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
111: M*/
115: PetscErrorCode SNESCreate_Test(SNES snes)
116: {
117: SNES_Test *neP;
121: snes->ops->setup = 0;
122: snes->ops->solve = SNESSolve_Test;
123: snes->ops->destroy = SNESDestroy_Test;
124: snes->ops->setfromoptions = SNESSetFromOptions_Test;
125: snes->ops->converged = 0;
126: snes->ops->view = 0;
128: ierr = PetscNew(SNES_Test,&neP);
129: PetscLogObjectMemory(snes,sizeof(SNES_Test));
130: snes->data = (void*)neP;
131: neP->complete_print = PETSC_FALSE;
132: return(0);
133: }