1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 rosenbrock1 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* Include "petsctao.h" so we can use TAO solvers. */ 4c4762a1bSJed Brown #include <petsctao.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown static char help[] = "This example demonstrates use of the TAO package to \n\ 7c4762a1bSJed Brown solve an unconstrained minimization problem on a single processor. We \n\ 8c4762a1bSJed Brown minimize the extended Rosenbrock function: \n\ 9c4762a1bSJed Brown sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\ 10c4762a1bSJed Brown or the chained Rosenbrock function:\n\ 11c4762a1bSJed Brown sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n"; 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* 14c4762a1bSJed Brown User-defined application context - contains data needed by the 15c4762a1bSJed Brown application-provided call-back routines that evaluate the function, 16c4762a1bSJed Brown gradient, and hessian. 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown typedef struct { 19c4762a1bSJed Brown PetscInt n; /* dimension */ 20c4762a1bSJed Brown PetscReal alpha; /* condition parameter */ 21c4762a1bSJed Brown PetscBool chained; 22c4762a1bSJed Brown } AppCtx; 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* -------------- User-defined routines ---------- */ 25c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 26c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 27c4762a1bSJed Brown 28c4762a1bSJed Brown int main(int argc,char **argv) 29c4762a1bSJed Brown { 30c4762a1bSJed Brown PetscReal zero=0.0; 31c4762a1bSJed Brown Vec x; /* solution vector */ 32c4762a1bSJed Brown Mat H; 33c4762a1bSJed Brown Tao tao; /* Tao solver context */ 34c4762a1bSJed Brown PetscBool flg, test_lmvm = PETSC_FALSE; 35c4762a1bSJed Brown PetscMPIInt size; /* number of processes running */ 36c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 37c4762a1bSJed Brown KSP ksp; 38c4762a1bSJed Brown PC pc; 39c4762a1bSJed Brown Mat M; 40c4762a1bSJed Brown Vec in, out, out2; 41c4762a1bSJed Brown PetscReal mult_solve_dist; 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Initialize TAO and PETSc */ 44*327415f7SBarry Smith PetscFunctionBeginUser; 459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 473c859ba3SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Initialize problem parameters */ 50c4762a1bSJed Brown user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE; 51c4762a1bSJed Brown /* Check for command line arguments to override defaults */ 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg)); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Allocate vectors for the solution and gradient */ 589566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.n,&x)); 599566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* The TAO code begins here */ 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Create TAO solver with desired solution method */ 649566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 659566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOLMVM)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Set solution vec and an initial guess */ 689566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero)); 699566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Set routines for function, gradient, hessian evaluation */ 729566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user)); 739566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,H,H,FormHessian,&user)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Test the LMVM matrix */ 76c4762a1bSJed Brown if (test_lmvm) { 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bqnktr")); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Check for TAO command line options */ 819566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 849566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Test the LMVM matrix */ 87c4762a1bSJed Brown if (test_lmvm) { 889566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 899566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 909566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(pc, &M)); 919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &in)); 929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out)); 939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out2)); 949566063dSJacob Faibussowitsch PetscCall(VecSet(in, 1.0)); 959566063dSJacob Faibussowitsch PetscCall(MatMult(M, in, out)); 969566063dSJacob Faibussowitsch PetscCall(MatSolve(M, out, out2)); 979566063dSJacob Faibussowitsch PetscCall(VecAXPY(out2, -1.0, in)); 989566063dSJacob Faibussowitsch PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist)); 99c4762a1bSJed Brown if (mult_solve_dist < 1.e-11) { 1009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-11\n")); 101c4762a1bSJed Brown } else if (mult_solve_dist < 1.e-6) { 1029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-6\n")); 103c4762a1bSJed Brown } else { 1049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 105c4762a1bSJed Brown } 1069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&in)); 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out2)); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 114c4762a1bSJed Brown 1159566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 116b122ec5aSJacob Faibussowitsch return 0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X). 122c4762a1bSJed Brown 123c4762a1bSJed Brown Input Parameters: 124c4762a1bSJed Brown . tao - the Tao context 125c4762a1bSJed Brown . X - input vector 126c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunctionGradient() 127c4762a1bSJed Brown 128c4762a1bSJed Brown Output Parameters: 129c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 130c4762a1bSJed Brown . f - function value 131c4762a1bSJed Brown 132c4762a1bSJed Brown Note: 133c4762a1bSJed Brown Some optimization methods ask for the function and the gradient evaluation 134c4762a1bSJed Brown at the same time. Evaluating both at once may be more efficient that 135c4762a1bSJed Brown evaluating each separately. 136c4762a1bSJed Brown */ 137c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr) 138c4762a1bSJed Brown { 139c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 140c4762a1bSJed Brown PetscInt i,nn=user->n/2; 141c4762a1bSJed Brown PetscReal ff=0,t1,t2,alpha=user->alpha; 142c4762a1bSJed Brown PetscScalar *g; 143c4762a1bSJed Brown const PetscScalar *x; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBeginUser; 146c4762a1bSJed Brown /* Get pointers to vector data */ 1479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetArray(G,&g)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* Compute G(X) */ 151c4762a1bSJed Brown if (user->chained) { 152c4762a1bSJed Brown g[0] = 0; 153c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 154c4762a1bSJed Brown t1 = x[i+1] - x[i]*x[i]; 155c4762a1bSJed Brown ff += PetscSqr(1 - x[i]) + alpha*t1*t1; 156c4762a1bSJed Brown g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]); 157c4762a1bSJed Brown g[i+1] = 2*alpha*t1; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } else { 160c4762a1bSJed Brown for (i=0; i<nn; i++) { 161c4762a1bSJed Brown t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i]; 162c4762a1bSJed Brown ff += alpha*t1*t1 + t2*t2; 163c4762a1bSJed Brown g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2; 164c4762a1bSJed Brown g[2*i+1] = 2*alpha*t1; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Restore vectors */ 1699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 1709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G,&g)); 171c4762a1bSJed Brown *f = ff; 172c4762a1bSJed Brown 1739566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(15.0*nn)); 174c4762a1bSJed Brown PetscFunctionReturn(0); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 180c4762a1bSJed Brown 181c4762a1bSJed Brown Input Parameters: 182c4762a1bSJed Brown . tao - the Tao context 183c4762a1bSJed Brown . x - input vector 184c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 185c4762a1bSJed Brown 186c4762a1bSJed Brown Output Parameters: 187c4762a1bSJed Brown . H - Hessian matrix 188c4762a1bSJed Brown 189c4762a1bSJed Brown Note: Providing the Hessian may not be necessary. Only some solvers 190c4762a1bSJed Brown require this matrix. 191c4762a1bSJed Brown */ 192c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 193c4762a1bSJed Brown { 194c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 195c4762a1bSJed Brown PetscInt i, ind[2]; 196c4762a1bSJed Brown PetscReal alpha=user->alpha; 197c4762a1bSJed Brown PetscReal v[2][2]; 198c4762a1bSJed Brown const PetscScalar *x; 199c4762a1bSJed Brown PetscBool assembled; 200c4762a1bSJed Brown 201c4762a1bSJed Brown PetscFunctionBeginUser; 202c4762a1bSJed Brown /* Zero existing matrix entries */ 2039566063dSJacob Faibussowitsch PetscCall(MatAssembled(H,&assembled)); 2049566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(H)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Get a pointer to vector data */ 2079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* Compute H(X) entries */ 210c4762a1bSJed Brown if (user->chained) { 2119566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 212c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 213c4762a1bSJed Brown PetscScalar t1 = x[i+1] - x[i]*x[i]; 214c4762a1bSJed Brown v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]); 215c4762a1bSJed Brown v[0][1] = 2*alpha*(-2*x[i]); 216c4762a1bSJed Brown v[1][0] = 2*alpha*(-2*x[i]); 217c4762a1bSJed Brown v[1][1] = 2*alpha*t1; 218c4762a1bSJed Brown ind[0] = i; ind[1] = i+1; 2199566063dSJacob Faibussowitsch PetscCall(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES)); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown } else { 222c4762a1bSJed Brown for (i=0; i<user->n/2; i++) { 223c4762a1bSJed Brown v[1][1] = 2*alpha; 224c4762a1bSJed Brown v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2; 225c4762a1bSJed Brown v[1][0] = v[0][1] = -4.0*alpha*x[2*i]; 226c4762a1bSJed Brown ind[0]=2*i; ind[1]=2*i+1; 2279566063dSJacob Faibussowitsch PetscCall(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES)); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown } 2309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* Assemble matrix */ 2339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 2349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 2359566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0*user->n/2.0)); 236c4762a1bSJed Brown PetscFunctionReturn(0); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /*TEST 240c4762a1bSJed Brown 241c4762a1bSJed Brown build: 242c4762a1bSJed Brown requires: !complex 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -tao_gatol 1.e-4 246c4762a1bSJed Brown requires: !single 247c4762a1bSJed Brown 248c4762a1bSJed Brown test: 249c4762a1bSJed Brown suffix: 2 250c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -tao_gatol 1.e-3 251c4762a1bSJed Brown 252c4762a1bSJed Brown test: 253c4762a1bSJed Brown suffix: 3 254c4762a1bSJed Brown args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4 255c4762a1bSJed Brown requires: !single 256c4762a1bSJed Brown 257c4762a1bSJed Brown test: 258c4762a1bSJed Brown suffix: 4 259c4762a1bSJed Brown args: -tao_smonitor -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4 260c4762a1bSJed Brown 261c4762a1bSJed Brown test: 262c4762a1bSJed Brown suffix: 5 263c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: 6 267c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -tao_gatol 1.e-4 268c4762a1bSJed Brown 269c4762a1bSJed Brown test: 270c4762a1bSJed Brown suffix: 7 271c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -tao_gatol 1.e-4 272c4762a1bSJed Brown 273c4762a1bSJed Brown test: 274c4762a1bSJed Brown suffix: 8 275c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 276c4762a1bSJed Brown 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown suffix: 9 279c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: 10 283c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 284c4762a1bSJed Brown 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown suffix: 11 287864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbroyden 288c4762a1bSJed Brown 289c4762a1bSJed Brown test: 290c4762a1bSJed Brown suffix: 12 291864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbroyden 292c4762a1bSJed Brown 293c4762a1bSJed Brown test: 294c4762a1bSJed Brown suffix: 13 295864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbroyden 296c4762a1bSJed Brown 297c4762a1bSJed Brown test: 298c4762a1bSJed Brown suffix: 14 299c4762a1bSJed Brown args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs 300c4762a1bSJed Brown 301c4762a1bSJed Brown test: 302c4762a1bSJed Brown suffix: 15 303c4762a1bSJed Brown args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: 16 307c4762a1bSJed Brown args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1 308c4762a1bSJed Brown 309c4762a1bSJed Brown test: 310c4762a1bSJed Brown suffix: 17 311c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnls 312c4762a1bSJed Brown 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: 18 315c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-4 -tao_type blmvm 316c4762a1bSJed Brown 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 19 319c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 320c4762a1bSJed Brown 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown suffix: 20 323c4762a1bSJed Brown args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor 324c4762a1bSJed Brown 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 21 327864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbroyden 328c4762a1bSJed Brown 32905579b36STristan Konolige test: 33005579b36STristan Konolige suffix: 22 33105579b36STristan Konolige args: -tao_max_it 1 -tao_converged_reason 33205579b36STristan Konolige 33305579b36STristan Konolige test: 33405579b36STristan Konolige suffix: 23 33505579b36STristan Konolige args: -tao_max_funcs 0 -tao_converged_reason 33605579b36STristan Konolige 33705579b36STristan Konolige test: 33805579b36STristan Konolige suffix: 24 33905579b36STristan Konolige args: -tao_gatol 10 -tao_converged_reason 34005579b36STristan Konolige 34105579b36STristan Konolige test: 34205579b36STristan Konolige suffix: 25 34305579b36STristan Konolige args: -tao_grtol 10 -tao_converged_reason 34405579b36STristan Konolige 34505579b36STristan Konolige test: 34605579b36STristan Konolige suffix: 26 34705579b36STristan Konolige args: -tao_gttol 10 -tao_converged_reason 34805579b36STristan Konolige 34905579b36STristan Konolige test: 35005579b36STristan Konolige suffix: 27 35105579b36STristan Konolige args: -tao_steptol 10 -tao_converged_reason 35205579b36STristan Konolige 35305579b36STristan Konolige test: 35405579b36STristan Konolige suffix: 28 35505579b36STristan Konolige args: -tao_fmin 10 -tao_converged_reason 35605579b36STristan Konolige 357c4762a1bSJed Brown TEST*/ 358