Actual source code: fdtest.c
1: #include "tao.h"
2: #ifdef TAO_USE_PETSC
3: #include "src/petsctao/application/taoapp/taoapp_petsc.h"
5: typedef struct {
6: PetscTruth check_gradient;
7: PetscTruth check_hessian;
8: PetscTruth complete_print;
9: } TAO_FD;
13: static int TaoSolve_FD(TAO_SOLVER tao, void *solver)
14: {
15: Mat A,Apre,B;
16: MatStructure flg;
17: Vec x,g1,g2;
18: PetscInt i;
19: PetscReal hcnorm,diffnorm;
20: int info;
21: TAO_FD *fd = (TAO_FD*)solver;
22: TaoApplication *tapp;
23: TaoPetscApplication *tpapp;
24: TAO_APPLICATION app;
25: MPI_Comm comm, gcomm;
26: TaoFunctionBegin;
27:
28: info = TaoGetApplication(tao,&tapp); CHKERRQ(info);
29: tpapp = dynamic_cast<TaoPetscApplication*>(tapp);
30: if (!tpapp) {
31: SETERRQ(1,"tao_fd_check solver only works with PETSc applications");
32: }
33: app = tpapp->papp;
34: info = tpapp->GetComm(&comm); CHKERRQ(info);
35: info = TaoAppGetSolutionVec(app,&x); CHKERRQ(info);
36: if (fd->check_gradient) {
37: info = PetscPrintf(comm,"Testing hand-coded gradient, if the ratio ||fd - hc|| / ||hc|| is\n"); CHKERRQ(info);
38: info = PetscPrintf(comm,"0 (1.e-8), the hand-coded gradient is probably correct.\n"); CHKERRQ(info);
40: if (!fd->complete_print) {
41: info = PetscPrintf(comm,"Run with -tao_test_display to show difference\n");CHKERRQ(info);
42: info = PetscPrintf(comm,"between hand-coded and finite difference gradient.\n");CHKERRQ(info);
43: }
44:
45: info = VecDuplicate(x,&g1); CHKERRQ(info);
46: info = VecDuplicate(x,&g2); CHKERRQ(info);
47: for (i=0; i<3; i++) {
48: if (i == 1) {info = VecSet(x,-1.0);CHKERRQ(info);}
49: else if (i == 2) {info = VecSet(x,1.0);CHKERRQ(info);}
50:
51: /* compute both versions of gradient */
52: info = TaoAppComputeGradient(app,x,g1); CHKERRQ(info);
53: info = TaoAppDefaultComputeGradient(app,x,g2,PETSC_NULL); CHKERRQ(info);
54: if (fd->complete_print) {
55: info = PetscPrintf(comm,"Finite difference gradient\n"); CHKERRQ(info);
56: info = PetscObjectGetComm((PetscObject)g2,&gcomm); CHKERRQ(info);
57: info = VecView(g2,PETSC_VIEWER_STDOUT_(gcomm)); CHKERRQ(info);
58: info = PetscPrintf(comm,"Hand-coded gradient\n");CHKERRQ(info);
59: info = PetscObjectGetComm((PetscObject)g1,&gcomm);CHKERRQ(info);
60: info = VecView(g1,PETSC_VIEWER_STDOUT_(gcomm));CHKERRQ(info);
61: info = PetscPrintf(comm,"\n");CHKERRQ(info);
62: }
63:
64: info = VecAXPY(g2,-1.0,g1); CHKERRQ(info); // g2 := g2-g1
65: info = VecNorm(g1,NORM_2,&hcnorm); CHKERRQ(info);
66: info = VecNorm(g2,NORM_2,&diffnorm); CHKERRQ(info);
68: info = PetscPrintf(comm,"ratio ||fd-hc||/||hc|| = %G, difference ||fd-hc|| = %G\n",diffnorm/hcnorm,diffnorm);CHKERRQ(info);
69: info = PetscPrintf(comm,"\n");CHKERRQ(info);
70: }
71: }
72: if (fd->check_hessian) {
73: info = TaoAppGetHessianMat(app,&A, &Apre); CHKERRQ(info);
74: if (A != Apre) {
75: SETERRQ(1,"Cannot test with alternative preconditioner");
76: }
77: info = PetscPrintf(comm,"Testing hand-coded hessian, if the ratio ||fd - hc|| / ||hc|| is\n"); CHKERRQ(info);
78: info = PetscPrintf(comm,"0 (1.e-8), the hand-coded hessian is probably correct.\n"); CHKERRQ(info);
80: if (!fd->complete_print) {
81: info = PetscPrintf(comm,"Run with -tao_test_display to show difference\n");CHKERRQ(info);
82: info = PetscPrintf(comm,"between hand-coded and finite difference hessian.\n");CHKERRQ(info);
83: }
84: for (i=0;i<3;i++) {
85: if (i==1) {
86: info = VecSet(x,-1.0); CHKERRQ(info);
87: } else if (i==2) {
88: info = VecSet(x,1.0); CHKERRQ(info);
89: }
90: info = TaoAppComputeHessian(app,x,&A,&A,&flg); CHKERRQ(info);
91: if (!i) {
92: info = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&B); CHKERRQ(info);
93: }
94: info = TaoAppDefaultComputeHessian(app, x, &B, &B, &flg, PETSC_NULL); CHKERRQ(info);
95: if (fd->complete_print) {
96: info = PetscObjectGetComm((PetscObject)A,&gcomm); CHKERRQ(info);
97: info = PetscPrintf(comm,"Finite difference hessian\n"); CHKERRQ(info);
98: info = MatView(B,PETSC_VIEWER_STDOUT_(gcomm)); CHKERRQ(info);
99: info = PetscPrintf(comm,"Hand-coded hessian\n"); CHKERRQ(info);
100: info = MatView(A,PETSC_VIEWER_STDOUT_(gcomm)); CHKERRQ(info);
101: info = PetscPrintf(comm,"\n"); CHKERRQ(info);
102: }
103: info = MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN); CHKERRQ(info);
104: info = MatNorm(B,NORM_FROBENIUS,&diffnorm); CHKERRQ(info);
105: info = MatNorm(A,NORM_FROBENIUS,&hcnorm); CHKERRQ(info);
107: info = PetscPrintf(comm,"ratio ||fd-hc||/||hc|| = %G, difference ||fd-hc|| = %G\n",diffnorm/hcnorm,diffnorm);CHKERRQ(info);
108: info = PetscPrintf(comm,"\n");CHKERRQ(info);
109: }
110: info = MatDestroy(B); CHKERRQ(info);
111: }
112:
114: TaoFunctionReturn(0);
115: // PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
117: }
121: static int TaoSetUp_FD(TAO_SOLVER tao, void *solver)
122: {
123: //TAO_FD *fdctx = (TAO_FD*)solver;
124: TaoVec *X, *G;
125: int info;
126: TaoFunctionBegin;
127: info = TaoGetSolution(tao,&X); CHKERRQ(info);
128: info = X->Clone(&G); CHKERRQ(info);
129: info = TaoSetLagrangianGradientVector(tao,G); CHKERRQ(info);
130: info = TaoCheckFG(tao); CHKERRQ(info);
131: TaoFunctionReturn(0);
132: }
136: static int TaoDestroy_FD(TAO_SOLVER tao, void *solver)
137: {
138: TaoFunctionBegin;
139: TaoFunctionReturn(0);
140: }
144: static int TaoSetOptions_FD(TAO_SOLVER tao, void *solver)
145: {
146: TAO_FD *fd = (TAO_FD*)solver;
147: int info;
148: TaoFunctionBegin;
149: info = PetscOptionsHead("Hand-coded gradient/hessian tester options");CHKERRQ(info);
150: info = PetscOptionsName("-tao_test_display","Display difference between approximate and handcoded hessian","None",&fd->complete_print);CHKERRQ(info);
151: info = PetscOptionsName("-tao_test_gradient","Test Hand-coded gradient against finite difference gradient","None",&fd->check_gradient);CHKERRQ(info);
152: info = PetscOptionsName("-tao_test_hessian","Test Hand-coded hessian against finite difference hessian","None",&fd->check_hessian);CHKERRQ(info);
153: if (fd->check_gradient == PETSC_FALSE && fd->check_hessian == PETSC_FALSE) {
154: fd->check_gradient = PETSC_TRUE;
155: }
156:
157: info = PetscOptionsTail();CHKERRQ(info);
159: TaoFunctionReturn(0);
160: }
164: int TaoCreate_FD(TAO_SOLVER tao)
165: {
166: TAO_FD *fd;
167: int info;
168: TaoFunctionBegin;
169: info = TaoNew(TAO_FD,&fd); CHKERRQ(info);
170: info = PetscLogObjectMemory(tao,sizeof(TAO_FD)); CHKERRQ(info);
171: info = TaoSetTaoSolveRoutine(tao,TaoSolve_FD, (void*)fd); CHKERRQ(info);
172: info = TaoSetTaoSetUpDownRoutines(tao,TaoSetUp_FD, TaoDestroy_FD); CHKERRQ(info);
173: info = TaoSetTaoOptionsRoutine(tao, TaoSetOptions_FD); CHKERRQ(info);
174: // info = TaoSetTaoViewRoutine(tao,TaoView_FD); CHKERRQ(info);
175: fd->check_gradient=PETSC_TRUE;
176: fd->check_hessian=PETSC_FALSE;
177: fd->complete_print=PETSC_FALSE;
178: TaoFunctionReturn(0);
179: }
183: #else
187: int TaoCreate_FD(TAO_SOLVER tao)
188: {
189: TaoFunctionBegin;
190: SETERRQ(1,"tao_fd_test solver requires TAO_APPLICATION");
191: }
194: #endif