Actual source code: snestest.c
1: #define PETSCSNES_DLL
3: #include src/snes/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: PetscScalar mone = -1.0,one = 1.0;
23: PetscReal nrm,gnorm;
24: SNES_Test *neP = (SNES_Test*)snes->data;
28: if (A != snes->jacobian_pre) {
29: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot test with alternative preconditioner");
30: }
32: PetscPrintf(snes->comm,"Testing hand-coded Jacobian, if the ratio is\n");
33: PetscPrintf(snes->comm,"O(1.e-8), the hand-coded Jacobian is probably correct.\n");
34: if (!neP->complete_print) {
35: PetscPrintf(snes->comm,"Run with -snes_test_display to show difference\n");
36: PetscPrintf(snes->comm,"of hand-coded and finite difference Jacobian.\n");
37: }
39: for (i=0; i<3; i++) {
40: if (i == 1) {VecSet(x,mone);}
41: else if (i == 2) {VecSet(x,one);}
42:
43: /* compute both versions of Jacobian */
44: SNESComputeJacobian(snes,x,&A,&A,&flg);
45: if (!i) {MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&B);}
46: SNESDefaultComputeJacobian(snes,x,&B,&B,&flg,snes->funP);
47: if (neP->complete_print) {
48: MPI_Comm comm;
49: PetscPrintf(snes->comm,"Finite difference Jacobian\n");
50: PetscObjectGetComm((PetscObject)B,&comm);
51: MatView(B,PETSC_VIEWER_STDOUT_(comm));
52: }
53: /* compare */
54: MatAXPY(B,mone,A,DIFFERENT_NONZERO_PATTERN);
55: MatNorm(B,NORM_FROBENIUS,&nrm);
56: MatNorm(A,NORM_FROBENIUS,&gnorm);
57: if (neP->complete_print) {
58: MPI_Comm comm;
59: PetscPrintf(snes->comm,"Hand-coded Jacobian\n");
60: PetscObjectGetComm((PetscObject)B,&comm);
61: MatView(A,PETSC_VIEWER_STDOUT_(comm));
62: }
63: PetscPrintf(snes->comm,"Norm of matrix ratio %g difference %g\n",nrm/gnorm,nrm);
64: }
65: MatDestroy(B);
66: /*
67: Return error code cause Jacobian not good
68: */
69: PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
70: }
71: /* ------------------------------------------------------------ */
74: PetscErrorCode SNESDestroy_Test(SNES snes)
75: {
77: return(0);
78: }
82: static PetscErrorCode SNESSetFromOptions_Test(SNES snes)
83: {
84: SNES_Test *ls = (SNES_Test *)snes->data;
89: PetscOptionsHead("Hand-coded Jacobian tester options");
90: PetscOptionsName("-snes_test_display","Display difference between approximate and handcoded Jacobian","None",&ls->complete_print);
91: PetscOptionsTail();
92: return(0);
93: }
95: /* ------------------------------------------------------------ */
99: PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_Test(SNES snes)
100: {
101: SNES_Test *neP;
105: snes->setup = 0;
106: snes->solve = SNESSolve_Test;
107: snes->destroy = SNESDestroy_Test;
108: snes->converged = SNESConverged_LS;
109: snes->setfromoptions = SNESSetFromOptions_Test;
111: ierr = PetscNew(SNES_Test,&neP);
112: PetscLogObjectMemory(snes,sizeof(SNES_Test));
113: snes->data = (void*)neP;
114: neP->complete_print = PETSC_FALSE;
115: return(0);
116: }