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: */
 15: int SNESSolve_Test(SNES snes)
 16: {
 17:   Mat          A = snes->jacobian,B;
 18:   Vec          x = snes->vec_sol;
 19:   int          ierr,i;
 20:   MatStructure flg;
 21:   PetscScalar  mone = -1.0,one = 1.0;
 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(&mone,x);}
 40:     else if (i == 2) {VecSet(&one,x);}
 41: 
 42:     /* compute both versions of Jacobian */
 43:     SNESComputeJacobian(snes,x,&A,&A,&flg);
 44:     if (!i) {MatConvert(A,MATSAME,&B);}
 45:     SNESDefaultComputeJacobian(snes,x,&B,&B,&flg,snes->funP);
 46:     if (neP->complete_print) {
 47:       MPI_Comm comm;
 48:       PetscPrintf(snes->comm,"Finite difference Jacobian\n");
 49:       PetscObjectGetComm((PetscObject)B,&comm);
 50:       MatView(B,PETSC_VIEWER_STDOUT_(comm));
 51:     }
 52:     /* compare */
 53:     MatAXPY(&mone,A,B,DIFFERENT_NONZERO_PATTERN);
 54:     MatNorm(B,NORM_FROBENIUS,&nrm);
 55:     MatNorm(A,NORM_FROBENIUS,&gnorm);
 56:     if (neP->complete_print) {
 57:       MPI_Comm comm;
 58:       PetscPrintf(snes->comm,"Hand-coded Jacobian\n");
 59:       PetscObjectGetComm((PetscObject)B,&comm);
 60:       MatView(A,PETSC_VIEWER_STDOUT_(comm));
 61:     }
 62:     PetscPrintf(snes->comm,"Norm of matrix ratio %g difference %g\n",nrm/gnorm,nrm);
 63:   }
 64:   MatDestroy(B);
 65:   /*
 66:          Return error code cause Jacobian not good
 67:   */
 68:   PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
 69: }
 70: /* ------------------------------------------------------------ */
 73: int SNESDestroy_Test(SNES snes)
 74: {
 76:   return(0);
 77: }

 81: static int SNESSetFromOptions_Test(SNES snes)
 82: {
 83:   SNES_Test  *ls = (SNES_Test *)snes->data;
 84:   int        ierr;


 88:   PetscOptionsHead("Hand-coded Jacobian tester options");
 89:     PetscOptionsName("-snes_test_display","Display difference between approximate and handcoded Jacobian","None",&ls->complete_print);
 90:   PetscOptionsTail();
 91:   return(0);
 92: }

 94: /* ------------------------------------------------------------ */
 95: EXTERN_C_BEGIN
 98: int SNESCreate_Test(SNES  snes)
 99: {
100:   SNES_Test *neP;

104:   snes->setup                = 0;
105:   snes->solve                = SNESSolve_Test;
106:   snes->destroy                = SNESDestroy_Test;
107:   snes->converged        = SNESConverged_LS;
108:   snes->setfromoptions  = SNESSetFromOptions_Test;

110:   ierr                        = PetscNew(SNES_Test,&neP);
111:   PetscLogObjectMemory(snes,sizeof(SNES_Test));
112:   snes->data            = (void*)neP;
113:   neP->complete_print   = PETSC_FALSE;
114:   return(0);
115: }
116: EXTERN_C_END