1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 2a9603a14SPatrick Farrell #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 3aaa7dc30SBarry Smith #include <../src/tao/bound/impls/blmvm/blmvm.h> 4a7e14dcfSSatish Balay 5cd929ea3SAlp Dener #define BLMVM_STEP_BFGS 0 6cd929ea3SAlp Dener #define BLMVM_STEP_GRAD 1 7cd929ea3SAlp Dener 8cd929ea3SAlp Dener #define BLMVM_AS_NONE 0 9cd929ea3SAlp Dener #define BLMVM_AS_BERTSEKAS 1 10cd929ea3SAlp Dener #define BLMVM_AS_SIZE 2 11cd929ea3SAlp Dener 12cd929ea3SAlp Dener static const char *BLMVM_AS_TYPE[64] = {"none", "bertsekas"}; 13cd929ea3SAlp Dener 14cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBLMVMEstimateActiveSet(Tao tao, PetscInt asType) 15cd929ea3SAlp Dener { 16cd929ea3SAlp Dener PetscErrorCode ierr; 17cd929ea3SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 18cd929ea3SAlp Dener 19cd929ea3SAlp Dener PetscFunctionBegin; 20cd929ea3SAlp Dener if (!tao->bounded) PetscFunctionReturn(0); 21cd929ea3SAlp Dener switch (asType) { 22cd929ea3SAlp Dener case BLMVM_AS_NONE: 23cd929ea3SAlp Dener ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr); 24cd929ea3SAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, blmP->unprojected_gradient, tao->XU, PETSC_TRUE, &blmP->inactive_idx);CHKERRQ(ierr); 25cd929ea3SAlp Dener ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr); 26cd929ea3SAlp Dener ierr = ISComplementVec(blmP->inactive_idx, tao->solution, &blmP->active_idx);CHKERRQ(ierr); 27cd929ea3SAlp Dener break; 28cd929ea3SAlp Dener 29cd929ea3SAlp Dener case BLMVM_AS_BERTSEKAS: 30cd929ea3SAlp Dener /* Use gradient descent to estimate the active set */ 31cd929ea3SAlp Dener ierr = VecCopy(blmP->unprojected_gradient, blmP->W);CHKERRQ(ierr); 32cd929ea3SAlp Dener ierr = VecScale(blmP->W, -1.0);CHKERRQ(ierr); 33cd929ea3SAlp Dener ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, blmP->unprojected_gradient, blmP->W, blmP->work, blmP->as_step, &blmP->as_tol, 34cd929ea3SAlp Dener &blmP->active_lower, &blmP->active_upper, &blmP->active_fixed, &blmP->active_idx, &blmP->inactive_idx);CHKERRQ(ierr); 35cd929ea3SAlp Dener break; 36cd929ea3SAlp Dener 37cd929ea3SAlp Dener default: 38cd929ea3SAlp Dener break; 39cd929ea3SAlp Dener } 40cd929ea3SAlp Dener PetscFunctionReturn(0); 41cd929ea3SAlp Dener } 42cd929ea3SAlp Dener 43cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBLMVMBoundStep(Tao tao, PetscInt asType, Vec step) 44cd929ea3SAlp Dener { 45cd929ea3SAlp Dener PetscErrorCode ierr; 46cd929ea3SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 47cd929ea3SAlp Dener 48cd929ea3SAlp Dener PetscFunctionBegin; 49cd929ea3SAlp Dener if (!tao->bounded) PetscFunctionReturn(0); 50cd929ea3SAlp Dener switch (asType) { 51cd929ea3SAlp Dener case BLMVM_AS_NONE: 52cd929ea3SAlp Dener ierr = VecISSet(step, blmP->active_idx, 0.0);CHKERRQ(ierr); 53cd929ea3SAlp Dener break; 54cd929ea3SAlp Dener 55cd929ea3SAlp Dener case BLMVM_AS_BERTSEKAS: 56cd929ea3SAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, blmP->active_lower, blmP->active_upper, blmP->active_fixed, 1.0, step);CHKERRQ(ierr); 57cd929ea3SAlp Dener break; 58cd929ea3SAlp Dener 59cd929ea3SAlp Dener default: 60cd929ea3SAlp Dener break; 61cd929ea3SAlp Dener } 62cd929ea3SAlp Dener PetscFunctionReturn(0); 63cd929ea3SAlp Dener } 64cd929ea3SAlp Dener 65a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 66cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSolve_BLMVM(Tao tao) 67a7e14dcfSSatish Balay { 68a7e14dcfSSatish Balay PetscErrorCode ierr; 69a7e14dcfSSatish Balay TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 70e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 71cd929ea3SAlp Dener PetscReal f, fold, gdx, gnorm, resnorm; 72*b2d8c577SAlp Dener PetscReal stepsize = 1.0; 73cd929ea3SAlp Dener PetscInt nDiff, nupdates, stepType = BLMVM_STEP_GRAD; 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay PetscFunctionBegin; 76a7e14dcfSSatish Balay /* Project initial point onto bounds */ 77a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 78cd929ea3SAlp Dener ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr); 79cd929ea3SAlp Dener if (tao->bounded) { 80a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 81cd929ea3SAlp Dener } 82a9603a14SPatrick Farrell 83a7e14dcfSSatish Balay /* Check convergence criteria */ 84a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr); 85a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 86a7e14dcfSSatish Balay 87cd929ea3SAlp Dener ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 883fcc4fd3STodd Munson if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 89a7e14dcfSSatish Balay 903ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 91cd929ea3SAlp Dener ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr); 92cd929ea3SAlp Dener ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr); 93cd929ea3SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 94cd929ea3SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr); 953ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 963ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 99cd929ea3SAlp Dener if (!blmP->recycle) { 100cd929ea3SAlp Dener blmP->bfgs = 0; 101a7e14dcfSSatish Balay blmP->grad = 0; 102cd929ea3SAlp Dener ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr); 103cd929ea3SAlp Dener } 104a7e14dcfSSatish Balay 105a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 1063ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 107cd929ea3SAlp Dener /* Estimate active set at the current iterate */ 108cd929ea3SAlp Dener ierr = TaoBLMVMEstimateActiveSet(tao, blmP->as_type);CHKERRQ(ierr); 109cd929ea3SAlp Dener 110a7e14dcfSSatish Balay /* Compute direction */ 111cd929ea3SAlp Dener if (blmP->H0) { 112cd929ea3SAlp Dener ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr); 113cd929ea3SAlp Dener stepType = BLMVM_STEP_BFGS; 114cd929ea3SAlp Dener } 115*b2d8c577SAlp Dener ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 116*b2d8c577SAlp Dener ierr = VecCopy(blmP->unprojected_gradient, blmP->red_grad);CHKERRQ(ierr); 117*b2d8c577SAlp Dener ierr = VecISSet(blmP->red_grad, blmP->active_idx, 0.0);CHKERRQ(ierr); 118*b2d8c577SAlp Dener ierr = MatSolve(blmP->M, blmP->red_grad, tao->stepdirection);CHKERRQ(ierr); 119*b2d8c577SAlp Dener ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 120*b2d8c577SAlp Dener ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr); 121cd929ea3SAlp Dener ierr = MatLMVMGetUpdateCount(blmP->M, &nupdates);CHKERRQ(ierr); 122cd929ea3SAlp Dener if (nupdates > 0) stepType = BLMVM_STEP_BFGS; 123a7e14dcfSSatish Balay 124a7e14dcfSSatish Balay /* Check for success (descent direction) */ 125*b2d8c577SAlp Dener ierr = VecDot(tao->stepdirection, blmP->red_grad, &gdx);CHKERRQ(ierr); 126*b2d8c577SAlp Dener if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) { 127a7e14dcfSSatish Balay /* Step is not descent or solve was not successful 128a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 129*b2d8c577SAlp Dener ierr = MatLMVMResetJ0(blmP->M);CHKERRQ(ierr); 130cd929ea3SAlp Dener ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr); 131a7e14dcfSSatish Balay ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 132cd929ea3SAlp Dener ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 133a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 134cd929ea3SAlp Dener ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr); 135*b2d8c577SAlp Dener stepType = BLMVM_STEP_GRAD; 136*b2d8c577SAlp Dener } 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay /* Perform the linesearch */ 139a7e14dcfSSatish Balay fold = f; 140a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr); 141a7e14dcfSSatish Balay ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr); 142a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 143a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 144a7e14dcfSSatish Balay 145cd929ea3SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && stepType != BLMVM_STEP_GRAD) { 146a7e14dcfSSatish Balay /* Linesearch failed 147a7e14dcfSSatish Balay Reset factors and use scaled (projected) gradient step */ 148a7e14dcfSSatish Balay f = fold; 149a7e14dcfSSatish Balay ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr); 150a7e14dcfSSatish Balay ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr); 151a7e14dcfSSatish Balay 152*b2d8c577SAlp Dener ierr = MatLMVMResetJ0(blmP->M);CHKERRQ(ierr); 153cd929ea3SAlp Dener ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr); 154a7e14dcfSSatish Balay ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 155cd929ea3SAlp Dener ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 156a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 157cd929ea3SAlp Dener ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr); 158cd929ea3SAlp Dener stepType = BLMVM_STEP_GRAD; 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 161a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 162cd929ea3SAlp Dener } 163a7e14dcfSSatish Balay 164a7e14dcfSSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 165cd929ea3SAlp Dener /* Line search failed on a gradient step, so just mark reason for divergence */ 166cd929ea3SAlp Dener f = fold; 167cd929ea3SAlp Dener ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr); 168cd929ea3SAlp Dener ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr); 169a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 170cd929ea3SAlp Dener } else { 171cd929ea3SAlp Dener /* LS found valid step, so tally the step type and compute projected gradient */ 172cd929ea3SAlp Dener switch (stepType) { 173cd929ea3SAlp Dener case BLMVM_STEP_BFGS: 174cd929ea3SAlp Dener ++blmP->bfgs; 175cd929ea3SAlp Dener break; 176cd929ea3SAlp Dener case BLMVM_STEP_GRAD: 177cd929ea3SAlp Dener ++blmP->grad; 178cd929ea3SAlp Dener break; 179cd929ea3SAlp Dener default: 180a7e14dcfSSatish Balay break; 181a7e14dcfSSatish Balay } 182cd929ea3SAlp Dener ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 183cd929ea3SAlp Dener ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 184cd929ea3SAlp Dener if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number"); 185a7e14dcfSSatish Balay } 186a7e14dcfSSatish Balay 187e4cb33bbSBarry Smith /* Check for converged */ 1888931d482SJason Sarich tao->niter++; 189cd929ea3SAlp Dener ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr); 190cd929ea3SAlp Dener ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr); 191cd929ea3SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 192cd929ea3SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr); 1933ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay PetscFunctionReturn(0); 196a7e14dcfSSatish Balay } 197a7e14dcfSSatish Balay 198cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSetup_BLMVM(Tao tao) 199a7e14dcfSSatish Balay { 200a7e14dcfSSatish Balay TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 201a7e14dcfSSatish Balay PetscInt n,N; 202a7e14dcfSSatish Balay PetscErrorCode ierr; 203a7e14dcfSSatish Balay 204a7e14dcfSSatish Balay PetscFunctionBegin; 205a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetup() */ 206a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr); 207a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr); 208302440fdSBarry Smith ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr); 209cd929ea3SAlp Dener ierr = VecDuplicate(tao->solution, &blmP->W);CHKERRQ(ierr); 210cd929ea3SAlp Dener ierr = VecDuplicate(tao->solution, &blmP->work);CHKERRQ(ierr); 211*b2d8c577SAlp Dener ierr = VecDuplicate(tao->solution, &blmP->red_grad);CHKERRQ(ierr); 212a7e14dcfSSatish Balay 213a7e14dcfSSatish Balay if (!tao->stepdirection) { 21453506e15SBarry Smith ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay if (!tao->gradient) { 217a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 220a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 221a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 222cd929ea3SAlp Dener ierr = MatSetSizes(blmP->M, n, n, N, N);CHKERRQ(ierr); 223cd929ea3SAlp Dener ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr); 224a9603a14SPatrick Farrell 225a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 226a9603a14SPatrick Farrell if (blmP->H0) { 227cd929ea3SAlp Dener ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr); 228a9603a14SPatrick Farrell } 229a7e14dcfSSatish Balay PetscFunctionReturn(0); 230a7e14dcfSSatish Balay } 231a7e14dcfSSatish Balay 232a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 233cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoDestroy_BLMVM(Tao tao) 234a7e14dcfSSatish Balay { 235a7e14dcfSSatish Balay TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 236a7e14dcfSSatish Balay PetscErrorCode ierr; 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay PetscFunctionBegin; 239a7e14dcfSSatish Balay if (tao->setupcalled) { 240a7e14dcfSSatish Balay ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr); 241a7e14dcfSSatish Balay ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr); 242a7e14dcfSSatish Balay ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr); 243cd929ea3SAlp Dener ierr = VecDestroy(&blmP->W);CHKERRQ(ierr); 244cd929ea3SAlp Dener ierr = VecDestroy(&blmP->work);CHKERRQ(ierr); 245*b2d8c577SAlp Dener ierr = VecDestroy(&blmP->red_grad);CHKERRQ(ierr); 246a7e14dcfSSatish Balay } 247cd929ea3SAlp Dener ierr = ISDestroy(&blmP->active_lower);CHKERRQ(ierr); 248cd929ea3SAlp Dener ierr = ISDestroy(&blmP->active_upper);CHKERRQ(ierr); 249cd929ea3SAlp Dener ierr = ISDestroy(&blmP->active_fixed);CHKERRQ(ierr); 250cd929ea3SAlp Dener ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr); 251cd929ea3SAlp Dener ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr); 252*b2d8c577SAlp Dener ierr = MatDestroy(&blmP->M);CHKERRQ(ierr); 253a9603a14SPatrick Farrell if (blmP->H0) { 254a9603a14SPatrick Farrell PetscObjectDereference((PetscObject)blmP->H0); 255a9603a14SPatrick Farrell } 256a9603a14SPatrick Farrell 257a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 258a7e14dcfSSatish Balay PetscFunctionReturn(0); 259a7e14dcfSSatish Balay } 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 262cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao) 263a7e14dcfSSatish Balay { 264cd929ea3SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 265a7e14dcfSSatish Balay PetscErrorCode ierr; 266a7e14dcfSSatish Balay 267a7e14dcfSSatish Balay PetscFunctionBegin; 2681a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr); 269cd929ea3SAlp Dener ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr); 270cd929ea3SAlp Dener ierr = PetscOptionsEList("-tao_blmvm_as_type","active set estimation method", "", BLMVM_AS_TYPE, BLMVM_AS_SIZE, BLMVM_AS_TYPE[blmP->as_type], &blmP->as_type,NULL);CHKERRQ(ierr); 271cd929ea3SAlp Dener ierr = PetscOptionsReal("-tao_blmvm_as_tol", "initial tolerance used when estimating actively bounded variables","",blmP->as_tol,&blmP->as_tol,NULL);CHKERRQ(ierr); 272cd929ea3SAlp Dener ierr = PetscOptionsReal("-tao_blmvm_as_step", "step length used when estimating actively bounded variables","",blmP->as_step,&blmP->as_step,NULL);CHKERRQ(ierr); 273a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 274*b2d8c577SAlp Dener ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr); 275a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 276a7e14dcfSSatish Balay PetscFunctionReturn(0); 277a7e14dcfSSatish Balay } 278a7e14dcfSSatish Balay 279a7e14dcfSSatish Balay 280a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 281cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 282a7e14dcfSSatish Balay { 283a7e14dcfSSatish Balay TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 284a7e14dcfSSatish Balay PetscBool isascii; 285a7e14dcfSSatish Balay PetscErrorCode ierr; 286cd929ea3SAlp Dener PetscInt recycled_its; 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay PetscFunctionBegin; 289a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 290a7e14dcfSSatish Balay if (isascii) { 291a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, " Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr); 292cd929ea3SAlp Dener if (lmP->recycle) { 293cd929ea3SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Recycle: on\n");CHKERRQ(ierr); 294cd929ea3SAlp Dener recycled_its = lmP->bfgs + lmP->grad; 295cd929ea3SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr); 296cd929ea3SAlp Dener } 297*b2d8c577SAlp Dener ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 298*b2d8c577SAlp Dener ierr = MatView(lmP->M, viewer);CHKERRQ(ierr); 299*b2d8c577SAlp Dener ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 300a7e14dcfSSatish Balay } 301a7e14dcfSSatish Balay PetscFunctionReturn(0); 302a7e14dcfSSatish Balay } 303a7e14dcfSSatish Balay 304cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 305a7e14dcfSSatish Balay { 306a7e14dcfSSatish Balay TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 307a7e14dcfSSatish Balay PetscErrorCode ierr; 308a7e14dcfSSatish Balay 309a7e14dcfSSatish Balay PetscFunctionBegin; 310441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 311a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 312a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 31353506e15SBarry Smith if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 314a7e14dcfSSatish Balay 315a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr); 316a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr); 317a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 318a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 319a7e14dcfSSatish Balay 320a7e14dcfSSatish Balay ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr); 321a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr); 322a7e14dcfSSatish Balay ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr); 323a7e14dcfSSatish Balay PetscFunctionReturn(0); 324a7e14dcfSSatish Balay } 325a7e14dcfSSatish Balay 326a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 3271522df2eSJason Sarich /*MC 3281522df2eSJason Sarich TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 3291522df2eSJason Sarich for nonlinear minimization with bound constraints. It is an extension 3301522df2eSJason Sarich of TAOLMVM 3311522df2eSJason Sarich 3321eb8069cSJason Sarich Level: beginner 3331522df2eSJason Sarich M*/ 334728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 335a7e14dcfSSatish Balay { 336a7e14dcfSSatish Balay TAO_BLMVM *blmP; 3378caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 338a7e14dcfSSatish Balay PetscErrorCode ierr; 339a7e14dcfSSatish Balay 340a7e14dcfSSatish Balay PetscFunctionBegin; 341a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BLMVM; 342a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BLMVM; 343a7e14dcfSSatish Balay tao->ops->view = TaoView_BLMVM; 344a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 345a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BLMVM; 346a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_BLMVM; 347a7e14dcfSSatish Balay 3483c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr); 349a9603a14SPatrick Farrell blmP->H0 = NULL; 350cd929ea3SAlp Dener blmP->as_step = 0.001; 351cd929ea3SAlp Dener blmP->as_tol = 0.001; 352cd929ea3SAlp Dener blmP->as_type = BLMVM_AS_BERTSEKAS; 353cd929ea3SAlp Dener blmP->recycle = PETSC_FALSE; 354cd929ea3SAlp Dener 355a7e14dcfSSatish Balay tao->data = (void*)blmP; 3566552cf8aSJason Sarich 3576552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3586552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 3596552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 360a7e14dcfSSatish Balay 361a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 36263b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 363a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 364441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 3655d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 366cd929ea3SAlp Dener 367cd929ea3SAlp Dener ierr = KSPInitializePackage();CHKERRQ(ierr); 368cd929ea3SAlp Dener ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr); 369cd929ea3SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr); 370cd929ea3SAlp Dener ierr = MatSetType(blmP->M, MATLBFGS);CHKERRQ(ierr); 371cd929ea3SAlp Dener ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr); 372a7e14dcfSSatish Balay PetscFunctionReturn(0); 373a7e14dcfSSatish Balay } 374a7e14dcfSSatish Balay 375cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMSetH0(Tao tao, Mat H0) 376a9603a14SPatrick Farrell { 377a9603a14SPatrick Farrell TAO_LMVM *lmP; 378a9603a14SPatrick Farrell TAO_BLMVM *blmP; 379b625d6c7SJed Brown TaoType type; 380a9603a14SPatrick Farrell PetscBool is_lmvm, is_blmvm; 381a9603a14SPatrick Farrell PetscErrorCode ierr; 382a9603a14SPatrick Farrell 383a9603a14SPatrick Farrell ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 384a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 385a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 386a9603a14SPatrick Farrell 387a9603a14SPatrick Farrell if (is_lmvm) { 388a9603a14SPatrick Farrell lmP = (TAO_LMVM *)tao->data; 389a9603a14SPatrick Farrell ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 390a9603a14SPatrick Farrell lmP->H0 = H0; 391a9603a14SPatrick Farrell } else if (is_blmvm) { 392a9603a14SPatrick Farrell blmP = (TAO_BLMVM *)tao->data; 393a9603a14SPatrick Farrell ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 394a9603a14SPatrick Farrell blmP->H0 = H0; 39574c66251SPatrick Farrell } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 396a9603a14SPatrick Farrell PetscFunctionReturn(0); 397a9603a14SPatrick Farrell } 398a9603a14SPatrick Farrell 399cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0(Tao tao, Mat *H0) 400a9603a14SPatrick Farrell { 401a9603a14SPatrick Farrell TAO_LMVM *lmP; 402a9603a14SPatrick Farrell TAO_BLMVM *blmP; 403b625d6c7SJed Brown TaoType type; 404a9603a14SPatrick Farrell PetscBool is_lmvm, is_blmvm; 405a9603a14SPatrick Farrell Mat M; 406a9603a14SPatrick Farrell 407a9603a14SPatrick Farrell PetscErrorCode ierr; 408a9603a14SPatrick Farrell 409a9603a14SPatrick Farrell ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 410a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 411a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 412a9603a14SPatrick Farrell 413a9603a14SPatrick Farrell if (is_lmvm) { 414a9603a14SPatrick Farrell lmP = (TAO_LMVM *)tao->data; 415a9603a14SPatrick Farrell M = lmP->M; 416a9603a14SPatrick Farrell } else if (is_blmvm) { 417a9603a14SPatrick Farrell blmP = (TAO_BLMVM *)tao->data; 418a9603a14SPatrick Farrell M = blmP->M; 41974c66251SPatrick Farrell } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 420cd929ea3SAlp Dener ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr); 421a9603a14SPatrick Farrell PetscFunctionReturn(0); 422a9603a14SPatrick Farrell } 423a9603a14SPatrick Farrell 424cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0KSP(Tao tao, KSP *ksp) 425a9603a14SPatrick Farrell { 426a9603a14SPatrick Farrell TAO_LMVM *lmP; 427a9603a14SPatrick Farrell TAO_BLMVM *blmP; 428b625d6c7SJed Brown TaoType type; 429a9603a14SPatrick Farrell PetscBool is_lmvm, is_blmvm; 430a9603a14SPatrick Farrell Mat M; 431a9603a14SPatrick Farrell PetscErrorCode ierr; 432a9603a14SPatrick Farrell 433a9603a14SPatrick Farrell ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 434a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 435a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 436a9603a14SPatrick Farrell 437a9603a14SPatrick Farrell if (is_lmvm) { 438a9603a14SPatrick Farrell lmP = (TAO_LMVM *)tao->data; 439a9603a14SPatrick Farrell M = lmP->M; 440a9603a14SPatrick Farrell } else if (is_blmvm) { 441a9603a14SPatrick Farrell blmP = (TAO_BLMVM *)tao->data; 442a9603a14SPatrick Farrell M = blmP->M; 44374c66251SPatrick Farrell } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 444cd929ea3SAlp Dener ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr); 445a9603a14SPatrick Farrell PetscFunctionReturn(0); 446a9603a14SPatrick Farrell } 447