Actual source code: fdtest.c

petsc-3.9.0 2018-04-07
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),\n");
 31:     PetscPrintf(comm,"if the ratio ||fd - hc|| / ||hc|| is\n");
 32:     PetscPrintf(comm,"O(1.e-8), the hand-coded gradient is probably correct.\n");

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

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

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

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

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

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

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

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

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

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

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

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

155:    Level: intermediate

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

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

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