Actual source code: rosenbrock1.c

  1: /* Program usage: mpiexec -n 1 rosenbrock1 [-help] [all TAO options] */

  3: /*  Include "petsctao.h" so we can use TAO solvers.  */
  4: #include <petsctao.h>
  5: #include "rosenbrock1.h" // defines AppCtx, AppCtxFormFunctionGradient(), and AppCtxFormHessian()

  7: static char help[] = "This example demonstrates use of the TAO package to \n\
  8: solve an unconstrained minimization problem on a single processor.  We \n\
  9: minimize the extended Rosenbrock function: \n\
 10:    sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\
 11: or the chained Rosenbrock function:\n\
 12:    sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n";

 14: /* -------------- User-defined routines ---------- */
 15: static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
 16: static PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);

 18: int main(int argc, char **argv)
 19: {
 20:   Vec         x; /* solution vector */
 21:   Mat         H, Hpre;
 22:   Tao         tao;  /* Tao solver context */
 23:   PetscMPIInt size; /* number of processes running */
 24:   AppCtx      user; /* user-defined application context */
 25:   MPI_Comm    comm;

 27:   /* Initialize TAO and PETSc */
 28:   PetscFunctionBeginUser;
 29:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
 30:   comm = PETSC_COMM_WORLD;
 31:   PetscCallMPI(MPI_Comm_size(comm, &size));
 32:   PetscCheck(size == 1, comm, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");

 34:   /* Initialize problem parameters */
 35:   PetscCall(AppCtxInitialize(comm, &user));

 37:   /* Allocate vector for the solution */
 38:   PetscCall(AppCtxCreateSolution(&user, &x));

 40:   /* Allocate the Hessian matrix */
 41:   PetscCall(AppCtxCreateHessianMatrices(&user, &H, &Hpre));

 43:   /* The TAO code begins here */

 45:   /* Create TAO solver with desired solution method */
 46:   PetscCall(TaoCreate(comm, &tao));
 47:   PetscCall(TaoSetType(tao, TAOLMVM));

 49:   /* Set solution vec and an initial guess */
 50:   PetscCall(VecZeroEntries(x));
 51:   PetscCall(TaoSetSolution(tao, x));

 53:   /* Set routines for function, gradient, hessian evaluation */
 54:   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user));
 55:   PetscCall(TaoSetHessian(tao, H, Hpre, FormHessian, &user));

 57:   /* Check for TAO command line options */
 58:   PetscCall(TaoSetFromOptions(tao));

 60:   /* SOLVE THE APPLICATION */
 61:   PetscCall(TaoSolve(tao));

 63:   /* Clean up */
 64:   PetscCall(AppCtxFinalize(&user, tao));
 65:   PetscCall(TaoDestroy(&tao));
 66:   PetscCall(VecDestroy(&x));
 67:   PetscCall(MatDestroy(&H));
 68:   PetscCall(MatDestroy(&Hpre));

 70:   PetscCall(PetscFinalize());
 71:   return 0;
 72: }

 74: /*
 75:   FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).

 77:   Input Parameters:
 78: + tao  - the Tao context
 79: . X    - input vector
 80: - ptr  - optional user-defined context, as set by TaoSetFunctionGradient()

 82:   Output Parameters:
 83: + f - function value
 84: - G - vector containing the newly evaluated gradient

 86:   Note:
 87:   Some optimization methods ask for the function and the gradient evaluation
 88:   at the same time.  Evaluating both at once may be more efficient that
 89:   evaluating each separately.
 90: */
 91: static PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
 92: {
 93:   AppCtx *user = (AppCtx *)ptr;

 95:   PetscFunctionBeginUser;
 96:   PetscCall(AppCtxFormFunctionGradient(user, X, f, G));
 97:   PetscFunctionReturn(PETSC_SUCCESS);
 98: }

100: /*
101:   FormHessian - Evaluates Hessian matrix.

103:   Input Parameters:
104: + tao   - the Tao context
105: . x     - input vector
106: - ptr   - optional user-defined context, as set by TaoSetHessian()

108:   Output Parameters:
109: + H     - Hessian matrix
110: - Hpre  - Preconditioning matrix

112:   Note:  Providing the Hessian may not be necessary.  Only some solvers
113:   require this matrix.
114: */
115: static PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
116: {
117:   AppCtx *user = (AppCtx *)ptr;

119:   PetscFunctionBeginUser;
120:   PetscCall(AppCtxFormHessian(user, X, H));
121:   // Manual Jacobi preconditioner for testing
122:   if (user->jacobi_pc) {
123:     Vec v;

125:     PetscCall(VecDuplicate(X, &v));
126:     PetscCall(MatGetDiagonal(H, v));
127:     PetscCall(MatZeroEntries(Hpre));
128:     PetscCall(MatDiagonalSet(Hpre, v, INSERT_VALUES));
129:     PetscCall(VecDestroy(&v));
130:   }
131:   PetscFunctionReturn(PETSC_SUCCESS);
132: }

134: /*TEST

136:    build:
137:      requires: !complex !single !quad !defined(PETSC_USE_64BIT_INDICES) !__float128

139:    test:
140:      args: -tao_monitor_short -tao_type nls -tao_gatol 1.e-4

142:    test:
143:      suffix: 2
144:      args: -tao_monitor_short -tao_type lmvm -tao_gatol 1.e-3

146:    test:
147:      suffix: 3
148:      args: -tao_monitor_short -tao_type ntr -tao_gatol 1.e-4

150:    test:
151:      suffix: 4
152:      args: -tao_monitor_short -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4

154:    test:
155:      suffix: 5
156:      args: -tao_monitor_short -tao_type bntr -tao_gatol 1.e-4

158:    test:
159:      suffix: 6
160:      args: -tao_monitor_short -tao_type bntl -tao_gatol 1.e-4

162:    test:
163:      suffix: 7
164:      args: -tao_monitor_short -tao_type bnls -tao_gatol 1.e-4

166:    test:
167:      suffix: 8
168:      args: -tao_monitor_short -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_bnk_cg_tao_monitor_short

170:    test:
171:      suffix: 9
172:      args: -tao_monitor_short -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_bnk_cg_tao_monitor_short

174:    test:
175:      suffix: 10
176:      args: -tao_monitor_short -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_bnk_cg_tao_monitor_short

178:    test:
179:      suffix: 11
180:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmbroyden

182:    test:
183:      suffix: 12
184:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbroyden

186:    test:
187:      suffix: 13
188:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbroyden -tao_bqnk_mat_lmvm_beta {{0.0 0.25 1.0}} -tao_bqnk_mat_lmvm_rho 0.75 -tao_bqnk_mat_lmvm_sigma_hist 2

190:    test:
191:      suffix: 14
192:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs -tao_bqnk_mat_lmvm_scale_type {{scalar diagonal}} -tao_bqnk_mat_lmvm_alpha {{0.0 0.25 0.5}} -tao_bqnk_mat_lmvm_theta 1.0

194:    test:
195:      suffix: 15
196:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp

198:    test:
199:      suffix: 16
200:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1

202:    test:
203:      suffix: 17
204:      args: -tao_monitor_short -tao_gatol 1e-4 -tao_type bqnls

206:    test:
207:      suffix: 18
208:      args: -tao_monitor_short -tao_gatol 1e-4 -tao_type blmvm

210:    test:
211:      suffix: 19
212:      args: -tao_monitor_short -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1

214:    test:
215:      suffix: 20
216:      args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor

218:    test:
219:      suffix: 21
220:      args: -test_lmvm -tao_type bqnktr -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbroyden

222:    test:
223:      suffix: 22
224:      args: -tao_max_it 1 -tao_converged_reason

226:    test:
227:      suffix: 23
228:      args: -tao_max_funcs 0 -tao_converged_reason

230:    test:
231:      suffix: 24
232:      args: -tao_gatol 10 -tao_converged_reason

234:    test:
235:      suffix: 25
236:      args: -tao_grtol 10 -tao_converged_reason

238:    test:
239:      suffix: 26
240:      args: -tao_gttol 10 -tao_converged_reason

242:    test:
243:      suffix: 27
244:      args: -tao_steptol 10 -tao_converged_reason

246:    test:
247:      suffix: 28
248:      args: -tao_fmin 10 -tao_converged_reason

250:    test:
251:      suffix: snes
252:      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -snes_atol 1.e-4 -pc_type none -tao_mf_hessian -ksp_type cg

254:    test:
255:      suffix: snes_ls_armijo
256:      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtonls -snes_atol 1.e-4 -pc_type none -tao_mf_hessian -snes_linesearch_monitor -snes_linesearch_order 1

258:    test:
259:      suffix: snes_tr_cgnegcurve_kmdc
260:      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -snes_atol 1.e-4 -pc_type none -ksp_type cg -snes_tr_kmdc 0.9 -ksp_converged_neg_curve -ksp_converged_reason

262:    test:
263:      suffix: snes_ls_lmvm
264:      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtonls -snes_atol 1.e-4 -pc_type lmvm -tao_mf_hessian

266:    test:
267:      suffix: add_terms_l2_no_pre
268:      args: -tao_type nls -tao_add_terms reg_ -reg_tao_term_type halfl2squared -tao_term_sum_reg_scale 0.3 -tao_monitor_short -tao_view ::ascii_info_detail

270:    test:
271:      suffix: add_terms_l1_no_pre
272:      args: -tao_type nls -tao_add_terms reg_ -reg_tao_term_type l1 -reg_tao_term_l1_epsilon 0.4 -tao_term_sum_reg_scale 0.3 -tao_monitor_short -tao_view ::ascii_info_detail

274:    test:
275:      suffix: hpre_is_not_h
276:      args: -tao_type nls -jacobi_pc 1 -tao_view ::ascii_info_detail -n 10

278:    test:
279:      suffix: param_none
280:      args: -tao_type nls -tao_add_terms reg_ -reg_tao_term_type halfl2squared -tao_term_sum_reg_scale 0.3 -tao_monitor_short
281:      args: -tao_view ::ascii_info_detail -reg_tao_term_parameters_mode none

283: TEST*/