Actual source code: snestest.c
1: /*$Id: snestest.c,v 1.58 2001/08/07 03:04:12 balay Exp $*/
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: */
13: int SNESSolve_Test(SNES snes,int *its)
14: {
15: Mat A = snes->jacobian,B;
16: Vec x = snes->vec_sol;
17: int ierr,i;
18: MatStructure flg;
19: PetscScalar mone = -1.0,one = 1.0;
20: PetscReal nrm,gnorm;
21: SNES_Test *neP = (SNES_Test*)snes->data;
24: *its = 0;
26: if (A != snes->jacobian_pre) {
27: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot test with alternative preconditioner");
28: }
30: PetscPrintf(snes->comm,"Testing hand-coded Jacobian, if the ratio isn");
31: PetscPrintf(snes->comm,"O(1.e-8), the hand-coded Jacobian is probably correct.n");
32: if (!neP->complete_print) {
33: PetscPrintf(snes->comm,"Run with -snes_test_display to show differencen");
34: PetscPrintf(snes->comm,"of hand-coded and finite difference Jacobian.n");
35: }
37: for (i=0; i<3; i++) {
38: if (i == 1) {VecSet(&mone,x);}
39: else if (i == 2) {VecSet(&one,x);}
40:
41: /* compute both versions of Jacobian */
42: SNESComputeJacobian(snes,x,&A,&A,&flg);
43: if (!i) {MatConvert(A,MATSAME,&B);}
44: SNESDefaultComputeJacobian(snes,x,&B,&B,&flg,snes->funP);
45: if (neP->complete_print) {
46: PetscPrintf(snes->comm,"Finite difference Jacobiann");
47: MatView(B,PETSC_VIEWER_STDOUT_WORLD);
48: }
49: /* compare */
50: MatAXPY(&mone,A,B,DIFFERENT_NONZERO_PATTERN);
51: MatNorm(B,NORM_FROBENIUS,&nrm);
52: MatNorm(A,NORM_FROBENIUS,&gnorm);
53: if (neP->complete_print) {
54: PetscPrintf(snes->comm,"Hand-coded Jacobiann");
55: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
56: }
57: PetscPrintf(snes->comm,"Norm of matrix ratio %g difference %gn",nrm/gnorm,nrm);
58: }
59: MatDestroy(B);
60: /*
61: Return error code cause Jacobian not good
62: */
63: PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
64: }
65: /* ------------------------------------------------------------ */
66: int SNESDestroy_Test(SNES snes)
67: {
69: return(0);
70: }
72: static int SNESSetFromOptions_Test(SNES snes)
73: {
74: SNES_Test *ls = (SNES_Test *)snes->data;
75: int ierr;
79: PetscOptionsHead("Hand-coded Jacobian tester options");
80: PetscOptionsName("-snes_test_display","Display difference between approximate and handcoded Jacobian","None",&ls->complete_print);
81: PetscOptionsTail();
82: return(0);
83: }
85: /* ------------------------------------------------------------ */
86: EXTERN_C_BEGIN
87: int SNESCreate_Test(SNES snes)
88: {
89: SNES_Test *neP;
93: snes->setup = 0;
94: snes->solve = SNESSolve_Test;
95: snes->destroy = SNESDestroy_Test;
96: snes->converged = SNESConverged_LS;
97: snes->setfromoptions = SNESSetFromOptions_Test;
99: ierr = PetscNew(SNES_Test,&neP);
100: PetscLogObjectMemory(snes,sizeof(SNES_Test));
101: snes->data = (void*)neP;
102: neP->complete_print = PETSC_FALSE;
103: return(0);
104: }
105: EXTERN_C_END