Actual source code: fdtest.c

petsc-3.8.3 2017-12-09
Report Typos and Errors
  1:  #include <petsc/private/taoimpl.h>

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

  9: /*
 10:      TaoSolve_Test - Tests whether a hand computed Hessian
 11:      matches one compute via finite differences.
 12: */
 13: PetscErrorCode TaoSolve_Test(Tao tao)
 14: {
 15:   Mat            A = tao->hessian,B;
 16:   Vec            x = tao->solution,g1,g2;
 18:   PetscInt       i;
 19:   PetscReal      nrm,gnorm,hcnorm,fdnorm,hcmax,fdmax,diffmax,diffnorm;
 20:   PetscScalar    dot;
 21:   MPI_Comm       comm;
 22:   Tao_Test        *fd = (Tao_Test*)tao->data;

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

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

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

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

 58:       VecNorm(g2,NORM_2,&fdnorm);
 59:       VecNorm(g1,NORM_2,&hcnorm);
 60:       VecNorm(g2,NORM_INFINITY,&fdmax);
 61:       VecNorm(g1,NORM_INFINITY,&hcmax);
 62:       VecDot(g1,g2,&dot);
 63:       VecAXPY(g2,-1.0,g1);
 64:       VecNorm(g2,NORM_2,&diffnorm);
 65:       VecNorm(g2,NORM_INFINITY,&diffmax);
 66:       PetscPrintf(comm,"||fd|| %g, ||hc|| = %g, angle cosine = (fd'hc)/||fd||||hc|| = %g\n", (double)fdnorm, (double)hcnorm, (double)(PetscRealPart(dot)/(fdnorm*hcnorm)));
 67:       PetscPrintf(comm,"2-norm ||fd-hc||/max(||hc||,||fd||) = %g, difference ||fd-hc|| = %g\n", (double)(diffnorm/PetscMax(hcnorm,fdnorm)), (double)diffnorm);
 68:       PetscPrintf(comm,"max-norm ||fd-hc||/max(||hc||,||fd||) = %g, difference ||fd-hc|| = %g\n", (double)(diffmax/PetscMax(hcmax,fdmax)), (double)diffmax);
 69:     }
 70:     VecDestroy(&g1);
 71:     VecDestroy(&g2);
 72:   }

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

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

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

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

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

130: static PetscErrorCode TaoSetFromOptions_Test(PetscOptionItems *PetscOptionsObject,Tao tao)
131: {
132:   Tao_Test        *fd = (Tao_Test *)tao->data;

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

147: /* ------------------------------------------------------------ */
148: /*C
149:       FD_TEST - Test hand-coded Hessian against finite difference Hessian

151:    Options Database:
152: .    -tao_test_display  Display difference between approximate and hand-coded Hessian

154:    Level: intermediate

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

158: */
159: PETSC_EXTERN PetscErrorCode  TaoCreate_Test(Tao  tao)
160: {
161:   Tao_Test        *fd;

165:   tao->ops->setup           = 0;
166:   tao->ops->solve           = TaoSolve_Test;
167:   tao->ops->destroy         = TaoDestroy_Test;
168:   tao->ops->setfromoptions  = TaoSetFromOptions_Test;
169:   tao->ops->view            = 0;
170:   PetscNewLog(tao,&fd);
171:   tao->data                 = (void*)fd;
172:   fd->complete_print        = PETSC_FALSE;
173:   fd->check_gradient        = PETSC_TRUE;
174:   fd->check_hessian         = PETSC_FALSE;
175:   return(0);
176: }