Actual source code: snestest.c

  2: #include <private/snesimpl.h>

  4: typedef struct {
  5:   PetscBool  complete_print;
  6: } SNES_Test;

  8: /*
  9:      SNESSolve_Test - Tests whether a hand computed Jacobian 
 10:      matches one compute via finite differences.
 11: */
 14: PetscErrorCode SNESSolve_Test(SNES snes)
 15: {
 16:   Mat            A = snes->jacobian,B;
 17:   Vec            x = snes->vec_sol,f = snes->vec_func;
 19:   PetscInt       i;
 20:   MatStructure   flg;
 21:   PetscReal      nrm,gnorm;
 22:   SNES_Test      *neP = (SNES_Test*)snes->data;


 26:   if (A != snes->jacobian_pre) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot test with alternative preconditioner");

 28:   PetscPrintf(((PetscObject)snes)->comm,"Testing hand-coded Jacobian, if the ratio is\n");
 29:   PetscPrintf(((PetscObject)snes)->comm,"O(1.e-8), the hand-coded Jacobian is probably correct.\n");
 30:   if (!neP->complete_print) {
 31:     PetscPrintf(((PetscObject)snes)->comm,"Run with -snes_test_display to show difference\n");
 32:     PetscPrintf(((PetscObject)snes)->comm,"of hand-coded and finite difference Jacobian.\n");
 33:   }

 35:   for (i=0; i<3; i++) {
 36:     static const char *const loc[] = {"user-defined state","constant state -1.0","constant state 1.0"};
 37:     if (i == 1) {VecSet(x,-1.0);}
 38:     else if (i == 2) {VecSet(x,1.0);}

 40:     /* evaluate the function at this point because SNESDefaultComputeJacobianColor() assumes that the function has been evaluated and put into snes->vec_func */
 41:     SNESComputeFunction(snes,x,f);
 42:     if (snes->domainerror) {
 43:       PetscPrintf(((PetscObject)snes)->comm,"Domain error at %s\n",loc[i]);
 44:       snes->domainerror = PETSC_FALSE;
 45:       continue;
 46:     }

 48:     /* compute both versions of Jacobian */
 49:     SNESComputeJacobian(snes,x,&A,&A,&flg);
 50:     if (!i) {MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&B);}
 51:     SNESDefaultComputeJacobian(snes,x,&B,&B,&flg,snes->funP);
 52:     if (neP->complete_print) {
 53:       MPI_Comm    comm;
 54:       PetscViewer viewer;
 55:       PetscPrintf(((PetscObject)snes)->comm,"Finite difference Jacobian (%s)\n",loc[i]);
 56:       PetscObjectGetComm((PetscObject)B,&comm);
 57:       PetscViewerASCIIGetStdout(comm,&viewer);
 58:       MatView(B,viewer);
 59:     }
 60:     /* compare */
 61:     MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
 62:     MatNorm(B,NORM_FROBENIUS,&nrm);
 63:     MatNorm(A,NORM_FROBENIUS,&gnorm);
 64:     if (neP->complete_print) {
 65:       MPI_Comm    comm;
 66:       PetscViewer viewer;
 67:       PetscPrintf(((PetscObject)snes)->comm,"Hand-coded Jacobian (%s)\n",loc[i]);
 68:       PetscObjectGetComm((PetscObject)B,&comm);
 69:       PetscViewerASCIIGetStdout(comm,&viewer);
 70:       MatView(A,viewer);
 71:       PetscPrintf(((PetscObject)snes)->comm,"Hand-coded minus finite difference Jacobian (%s)\n",loc[i]);
 72:       MatView(B,viewer);
 73:     }
 74:     if (!gnorm) gnorm = 1; /* just in case */
 75:     PetscPrintf(((PetscObject)snes)->comm,"Norm of matrix ratio %G difference %G (%s)\n",nrm/gnorm,nrm,loc[i]);
 76:   }
 77:   MatDestroy(&B);
 78:   /*
 79:          Return error code cause Jacobian not good
 80:   */
 81:   PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
 82: }
 83: /* ------------------------------------------------------------ */
 86: PetscErrorCode SNESDestroy_Test(SNES snes)
 87: {
 90:   PetscFree(snes->data);
 91:   return(0);
 92: }

 96: static PetscErrorCode SNESSetFromOptions_Test(SNES snes)
 97: {
 98:   SNES_Test      *ls = (SNES_Test *)snes->data;


103:   PetscOptionsHead("Hand-coded Jacobian tester options");
104:   PetscOptionsBool("-snes_test_display","Display difference between hand-coded and finite difference Jacobians","None",ls->complete_print,&ls->complete_print,PETSC_NULL);
105:   PetscOptionsTail();
106:   return(0);
107: }

109: /* ------------------------------------------------------------ */
110: /*MC
111:       SNESTEST - Test hand-coded Jacobian against finite difference Jacobian

113:    Options Database:
114: .    -snes_test_display  Display difference between approximate and hand-coded Jacobian

116:    Level: intermediate

118: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR

120: M*/
124: PetscErrorCode  SNESCreate_Test(SNES  snes)
125: {
126:   SNES_Test      *neP;

130:   snes->ops->solve           = SNESSolve_Test;
131:   snes->ops->destroy         = SNESDestroy_Test;
132:   snes->ops->setfromoptions  = SNESSetFromOptions_Test;
133:   snes->ops->view            = 0;
134:   snes->ops->setup           = 0;
135:   snes->ops->reset           = 0;

137:   PetscNewLog(snes,SNES_Test,&neP);
138:   snes->data            = (void*)neP;
139:   neP->complete_print   = PETSC_FALSE;
140:   return(0);
141: }