Actual source code: fdtest.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <petsc-private/taoimpl.h> /*I "petsctao.h" I*/

  3: typedef struct {
  4:   PetscBool  check_gradient;
  5:   PetscBool  check_hessian;
  6:   PetscBool  complete_print;
  7: } FD_Test;

  9: /*
 10:      TaoSolve_FD - Tests whether a hand computed Hessian
 11:      matches one compute via finite differences.
 12: */
 15: PetscErrorCode TaoSolve_FD(Tao tao)
 16: {
 17:   Mat            A = tao->hessian,B;
 18:   Vec            x = tao->solution,g1,g2;
 20:   PetscInt       i;
 21:   MatStructure   flg;
 22:   PetscReal      nrm,gnorm,hcnorm,fdnorm;
 23:   MPI_Comm       comm;
 24:   FD_Test        *fd = (FD_Test*)tao->data;

 27:   comm = ((PetscObject)tao)->comm;
 28:   if (fd->check_gradient) {
 29:     VecDuplicate(x,&g1);
 30:     VecDuplicate(x,&g2);

 32:     PetscPrintf(comm,"Testing hand-coded gradient (hc) against finite difference gradient (fd), if the ratio ||fd - hc|| / ||hc|| is\n");
 33:     PetscPrintf(comm,"0 (1.e-8), the hand-coded gradient is probably correct.\n");

 35:     if (!fd->complete_print) {
 36:       PetscPrintf(comm,"Run with -tao_fd_test_display to show difference\n");
 37:       PetscPrintf(comm,"between hand-coded and finite difference gradient.\n");
 38:     }
 39:     for (i=0; i<3; i++) {
 40:       if (i == 1) {VecSet(x,-1.0);}
 41:       else if (i == 2) {VecSet(x,1.0);}

 43:       /* Compute both version of gradient */
 44:       TaoComputeGradient(tao,x,g1);
 45:       TaoDefaultComputeGradient(tao,x,g2,NULL);
 46:       if (fd->complete_print) {
 47:         MPI_Comm gcomm;
 48:         PetscViewer viewer;
 49:         PetscPrintf(comm,"Finite difference gradient\n");
 50:         PetscObjectGetComm((PetscObject)g2,&gcomm);
 51:         PetscViewerASCIIGetStdout(gcomm,&viewer);
 52:         VecView(g2,viewer);
 53:         PetscPrintf(comm,"Hand-coded gradient\n");
 54:         PetscObjectGetComm((PetscObject)g1,&gcomm);
 55:         PetscViewerASCIIGetStdout(gcomm,&viewer);
 56:         VecView(g1,viewer);
 57:         PetscPrintf(comm,"\n");
 58:       }

 60:       VecAXPY(g2,-1.0,g1);
 61:       VecNorm(g1,NORM_2,&hcnorm);
 62:       VecNorm(g2,NORM_2,&fdnorm);

 64:       if (!hcnorm) hcnorm=1.0e-20;
 65:       PetscPrintf(comm,"ratio ||fd-hc||/||hc|| = %g, difference ||fd-hc|| = %g\n", (double)(fdnorm/hcnorm), (double)fdnorm);

 67:     }
 68:     VecDestroy(&g1);
 69:     VecDestroy(&g2);
 70:   }

 72:   if (fd->check_hessian) {
 73:     if (A != tao->hessian_pre) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot test with alternative preconditioner");

 75:     PetscPrintf(comm,"Testing hand-coded Hessian (hc) against finite difference Hessian (fd). If the ratio is\n");
 76:     PetscPrintf(comm,"O (1.e-8), the hand-coded Hessian is probably correct.\n");

 78:     if (!fd->complete_print) {
 79:       PetscPrintf(comm,"Run with -tao_fd_test_display to show difference\n");
 80:       PetscPrintf(comm,"of hand-coded and finite difference Hessian.\n");
 81:     }
 82:     for (i=0;i<3;i++) {
 83:       /* compute both versions of Hessian */
 84:       TaoComputeHessian(tao,x,&A,&A,&flg);
 85:       if (!i) {MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&B);}
 86:       TaoDefaultComputeHessian(tao,x,&B,&B,&flg,tao->user_hessP);
 87:       if (fd->complete_print) {
 88:         MPI_Comm    bcomm;
 89:         PetscViewer viewer;
 90:         PetscPrintf(comm,"Finite difference Hessian\n");
 91:         PetscObjectGetComm((PetscObject)B,&bcomm);
 92:         PetscViewerASCIIGetStdout(bcomm,&viewer);
 93:         MatView(B,viewer);
 94:       }
 95:       /* compare */
 96:       MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
 97:       MatNorm(B,NORM_FROBENIUS,&nrm);
 98:       MatNorm(A,NORM_FROBENIUS,&gnorm);
 99:       if (fd->complete_print) {
100:         MPI_Comm    hcomm;
101:         PetscViewer viewer;
102:         PetscPrintf(comm,"Hand-coded Hessian\n");
103:         PetscObjectGetComm((PetscObject)B,&hcomm);
104:         PetscViewerASCIIGetStdout(hcomm,&viewer);
105:         MatView(A,viewer);
106:         PetscPrintf(comm,"Hand-coded minus finite difference Hessian\n");
107:         MatView(B,viewer);
108:       }
109:       if (!gnorm) gnorm = 1.0e-20;
110:       PetscPrintf(comm,"ratio ||fd-hc||/||hc|| = %g, difference ||fd-hc|| = %g\n",(double)(nrm/gnorm),(double)nrm);
111:     }

113:     MatDestroy(&B);
114:   }
115:   tao->reason = TAO_CONVERGED_USER;
116:   return(0);
117: }
118: /* ------------------------------------------------------------ */
121: PetscErrorCode TaoDestroy_FD(Tao tao)
122: {

126:   PetscFree(tao->data);
127:   return(0);
128: }

132: static PetscErrorCode TaoSetFromOptions_FD(Tao tao)
133: {
134:   FD_Test        *fd = (FD_Test *)tao->data;

138:   PetscOptionsHead("Hand-coded Hessian tester options");
139:   PetscOptionsBool("-tao_fd_test_display","Display difference between hand-coded and finite difference Hessians","None",fd->complete_print,&fd->complete_print,NULL);
140:   PetscOptionsBool("-tao_fd_test_gradient","Test Hand-coded gradient against finite-difference gradient","None",fd->check_gradient,&fd->check_gradient,NULL);
141:   PetscOptionsBool("-tao_fd_test_hessian","Test Hand-coded hessian against finite-difference hessian","None",fd->check_hessian,&fd->check_hessian,NULL);
142:   if (fd->check_gradient == PETSC_FALSE && fd->check_hessian == PETSC_FALSE) {
143:     fd->check_gradient = PETSC_TRUE;
144:   }
145:   PetscOptionsTail();
146:   return(0);
147: }

149: /* ------------------------------------------------------------ */
150: /*C
151:       FD_TEST - Test hand-coded Hessian against finite difference Hessian

153:    Options Database:
154: .    -tao_fd_test_display  Display difference between approximate and hand-coded Hessian

156:    Level: intermediate

158: .seealso:  TaoCreate(), TaoSetType()

160: */
161: EXTERN_C_BEGIN
164: PetscErrorCode  TaoCreate_FD(Tao  tao)
165: {
166:   FD_Test        *fd;

170:   tao->ops->setup            = 0;
171:   tao->ops->solve            = TaoSolve_FD;
172:   tao->ops->destroy          = TaoDestroy_FD;
173:   tao->ops->setfromoptions  = TaoSetFromOptions_FD;
174:   tao->ops->view            = 0;
175:   PetscNewLog(tao,&fd);
176:   tao->data     = (void*)fd;
177:   fd->complete_print   = PETSC_FALSE;
178:   fd->check_gradient = PETSC_TRUE;
179:   fd->check_hessian = PETSC_FALSE;
180:   return(0);
181: }
182: EXTERN_C_END