1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h> 3a9603a14SPatrick Farrell #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 4aaa7dc30SBarry Smith #include <../src/tao/bound/impls/blmvm/blmvm.h> 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 7a7e14dcfSSatish Balay #undef __FUNCT__ 8a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_BLMVM" 9441846f8SBarry Smith static PetscErrorCode TaoSolve_BLMVM(Tao tao) 10a7e14dcfSSatish Balay { 11a7e14dcfSSatish Balay PetscErrorCode ierr; 12a7e14dcfSSatish Balay TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 13e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 14e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 15a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 16a7e14dcfSSatish Balay PetscReal stepsize = 1.0,delta; 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay PetscFunctionBegin; 19a7e14dcfSSatish Balay /* Project initial point onto bounds */ 20a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 21a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 22a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 23a7e14dcfSSatish Balay 24a9603a14SPatrick Farrell 25a7e14dcfSSatish Balay /* Check convergence criteria */ 26a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr); 27a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 28a7e14dcfSSatish Balay 29a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 303fcc4fd3STodd Munson if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 31a7e14dcfSSatish Balay 328931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr); 3353506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 34a7e14dcfSSatish Balay 35a7e14dcfSSatish Balay /* Set initial scaling for the function */ 36a7e14dcfSSatish Balay if (f != 0.0) { 37a7e14dcfSSatish Balay delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm); 3853506e15SBarry Smith } else { 39a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 40a7e14dcfSSatish Balay } 41a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr); 42*1feb5a2fSTodd Munson ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr); 43a7e14dcfSSatish Balay 44a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 45a7e14dcfSSatish Balay blmP->grad = 0; 46a7e14dcfSSatish Balay blmP->reset = 0; 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 49a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 50a7e14dcfSSatish Balay /* Compute direction */ 51a7e14dcfSSatish Balay ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 52a7e14dcfSSatish Balay ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 53a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 54a7e14dcfSSatish Balay 55a7e14dcfSSatish Balay /* Check for success (descent direction) */ 56a7e14dcfSSatish Balay ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr); 57a7e14dcfSSatish Balay if (gdx <= 0) { 58a7e14dcfSSatish Balay /* Step is not descent or solve was not successful 59a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 60a7e14dcfSSatish Balay ++blmP->grad; 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay if (f != 0.0) { 63a7e14dcfSSatish Balay delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm); 6453506e15SBarry Smith } else { 65a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 66a7e14dcfSSatish Balay } 67a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr); 68a7e14dcfSSatish Balay ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr); 69a7e14dcfSSatish Balay ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 70a7e14dcfSSatish Balay ierr = MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 71a7e14dcfSSatish Balay } 72a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 73a7e14dcfSSatish Balay 74a7e14dcfSSatish Balay /* Perform the linesearch */ 75a7e14dcfSSatish Balay fold = f; 76a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr); 77a7e14dcfSSatish Balay ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 79a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 80a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 81a7e14dcfSSatish Balay 82a7e14dcfSSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 83a7e14dcfSSatish Balay /* Linesearch failed 84a7e14dcfSSatish Balay Reset factors and use scaled (projected) gradient step */ 85a7e14dcfSSatish Balay ++blmP->reset; 86a7e14dcfSSatish Balay 87a7e14dcfSSatish Balay f = fold; 88a7e14dcfSSatish Balay ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr); 89a7e14dcfSSatish Balay ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr); 90a7e14dcfSSatish Balay 91a7e14dcfSSatish Balay if (f != 0.0) { 92a7e14dcfSSatish Balay delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm); 9353506e15SBarry Smith } else { 94a7e14dcfSSatish Balay delta = 2.0/ (gnorm*gnorm); 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr); 97a7e14dcfSSatish Balay ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr); 98a7e14dcfSSatish Balay ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr); 99a7e14dcfSSatish Balay ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 100a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 101a7e14dcfSSatish Balay 1023fcc4fd3STodd Munson /* This may be incorrect; linesearch has values for stepmax and stepmin 103a7e14dcfSSatish Balay that should be reset. */ 104302440fdSBarry Smith ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 105a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr); 106a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 107a7e14dcfSSatish Balay 108a7e14dcfSSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 109a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 110a7e14dcfSSatish Balay break; 111a7e14dcfSSatish Balay } 112a7e14dcfSSatish Balay } 113a7e14dcfSSatish Balay 114e4cb33bbSBarry Smith /* Check for converged */ 115a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 116a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 117a7e14dcfSSatish Balay 118a7e14dcfSSatish Balay 11953506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number"); 1208931d482SJason Sarich tao->niter++; 1218931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr); 122a7e14dcfSSatish Balay } 123a7e14dcfSSatish Balay PetscFunctionReturn(0); 124a7e14dcfSSatish Balay } 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay #undef __FUNCT__ 127a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_BLMVM" 128441846f8SBarry Smith static PetscErrorCode TaoSetup_BLMVM(Tao tao) 129a7e14dcfSSatish Balay { 130a7e14dcfSSatish Balay TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 131a7e14dcfSSatish Balay PetscInt n,N; 132a7e14dcfSSatish Balay PetscErrorCode ierr; 133a9603a14SPatrick Farrell KSP H0ksp; 134a7e14dcfSSatish Balay 135a7e14dcfSSatish Balay PetscFunctionBegin; 136a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetup() */ 137a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr); 138a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr); 139302440fdSBarry Smith ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr); 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay if (!tao->stepdirection) { 14253506e15SBarry Smith ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 143a7e14dcfSSatish Balay } 144a7e14dcfSSatish Balay if (!tao->gradient) { 145a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 146a7e14dcfSSatish Balay } 147a7e14dcfSSatish Balay if (!tao->XL) { 148a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr); 149e270355aSBarry Smith ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr); 150a7e14dcfSSatish Balay } 151a7e14dcfSSatish Balay if (!tao->XU) { 152a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr); 153e270355aSBarry Smith ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr); 154a7e14dcfSSatish Balay } 155a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 156a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 157a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 158a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);CHKERRQ(ierr); 159a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(blmP->M,tao->solution);CHKERRQ(ierr); 160a9603a14SPatrick Farrell 161a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 162a9603a14SPatrick Farrell if (blmP->H0) { 163a9603a14SPatrick Farrell const char *prefix; 164a9603a14SPatrick Farrell PC H0pc; 165a9603a14SPatrick Farrell 166a9603a14SPatrick Farrell ierr = MatLMVMSetH0(blmP->M, blmP->H0);CHKERRQ(ierr); 167a9603a14SPatrick Farrell ierr = MatLMVMGetH0KSP(blmP->M, &H0ksp);CHKERRQ(ierr); 168a9603a14SPatrick Farrell 169a9603a14SPatrick Farrell ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr); 170a9603a14SPatrick Farrell ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr); 171a9603a14SPatrick Farrell ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr); 172a9603a14SPatrick Farrell ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr); 173a9603a14SPatrick Farrell ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc, "tao_h0_");CHKERRQ(ierr); 174a9603a14SPatrick Farrell 175a9603a14SPatrick Farrell ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr); 176a9603a14SPatrick Farrell ierr = KSPSetUp(H0ksp);CHKERRQ(ierr); 177a9603a14SPatrick Farrell } 178a9603a14SPatrick Farrell 179a7e14dcfSSatish Balay PetscFunctionReturn(0); 180a7e14dcfSSatish Balay } 181a7e14dcfSSatish Balay 182a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 183a7e14dcfSSatish Balay #undef __FUNCT__ 184a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_BLMVM" 185441846f8SBarry Smith static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 186a7e14dcfSSatish Balay { 187a7e14dcfSSatish Balay TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 188a7e14dcfSSatish Balay PetscErrorCode ierr; 189a7e14dcfSSatish Balay 190a7e14dcfSSatish Balay PetscFunctionBegin; 191a7e14dcfSSatish Balay if (tao->setupcalled) { 192a7e14dcfSSatish Balay ierr = MatDestroy(&blmP->M);CHKERRQ(ierr); 193a7e14dcfSSatish Balay ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr); 194a7e14dcfSSatish Balay ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr); 195a7e14dcfSSatish Balay ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr); 196a7e14dcfSSatish Balay } 197a9603a14SPatrick Farrell 198a9603a14SPatrick Farrell if (blmP->H0) { 199a9603a14SPatrick Farrell PetscObjectDereference((PetscObject)blmP->H0); 200a9603a14SPatrick Farrell } 201a9603a14SPatrick Farrell 202a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 203a7e14dcfSSatish Balay PetscFunctionReturn(0); 204a7e14dcfSSatish Balay } 205a7e14dcfSSatish Balay 206a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 207a7e14dcfSSatish Balay #undef __FUNCT__ 208a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_BLMVM" 2094416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao) 210a7e14dcfSSatish Balay { 211a7e14dcfSSatish Balay PetscErrorCode ierr; 212a7e14dcfSSatish Balay 213a7e14dcfSSatish Balay PetscFunctionBegin; 2141a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr); 215a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 216a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 217a7e14dcfSSatish Balay PetscFunctionReturn(0); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay 221a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 222a7e14dcfSSatish Balay #undef __FUNCT__ 223a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_BLMVM" 224441846f8SBarry Smith static int TaoView_BLMVM(Tao tao, PetscViewer viewer) 225a7e14dcfSSatish Balay { 226a7e14dcfSSatish Balay TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 227a7e14dcfSSatish Balay PetscBool isascii; 228a7e14dcfSSatish Balay PetscErrorCode ierr; 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay PetscFunctionBegin; 231a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 232a7e14dcfSSatish Balay if (isascii) { 233a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 234a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr); 235a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 236a7e14dcfSSatish Balay } 237a7e14dcfSSatish Balay PetscFunctionReturn(0); 238a7e14dcfSSatish Balay } 239a7e14dcfSSatish Balay 240a7e14dcfSSatish Balay #undef __FUNCT__ 241a7e14dcfSSatish Balay #define __FUNCT__ "TaoComputeDual_BLMVM" 242441846f8SBarry Smith static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 243a7e14dcfSSatish Balay { 244a7e14dcfSSatish Balay TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 245a7e14dcfSSatish Balay PetscErrorCode ierr; 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay PetscFunctionBegin; 248441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 249a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 250a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 25153506e15SBarry 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"); 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr); 254a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr); 255a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 256a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 257a7e14dcfSSatish Balay 258a7e14dcfSSatish Balay ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr); 259a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr); 260a7e14dcfSSatish Balay ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr); 261a7e14dcfSSatish Balay PetscFunctionReturn(0); 262a7e14dcfSSatish Balay } 263a7e14dcfSSatish Balay 264a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 2651522df2eSJason Sarich /*MC 2661522df2eSJason Sarich TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 2671522df2eSJason Sarich for nonlinear minimization with bound constraints. It is an extension 2681522df2eSJason Sarich of TAOLMVM 2691522df2eSJason Sarich 2701522df2eSJason Sarich Options Database Keys: 2711522df2eSJason Sarich + -tao_lmm_vectors - number of vectors to use for approximation 2721522df2eSJason Sarich . -tao_lmm_scale_type - "none","scalar","broyden" 2731522df2eSJason Sarich . -tao_lmm_limit_type - "none","average","relative","absolute" 2741522df2eSJason Sarich . -tao_lmm_rescale_type - "none","scalar","gl" 2751522df2eSJason Sarich . -tao_lmm_limit_mu - mu limiting factor 2761522df2eSJason Sarich . -tao_lmm_limit_nu - nu limiting factor 2771522df2eSJason Sarich . -tao_lmm_delta_min - minimum delta value 2781522df2eSJason Sarich . -tao_lmm_delta_max - maximum delta value 2791522df2eSJason Sarich . -tao_lmm_broyden_phi - phi factor for Broyden scaling 2801522df2eSJason Sarich . -tao_lmm_scalar_alpha - alpha factor for scalar scaling 2811522df2eSJason Sarich . -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal 2821522df2eSJason Sarich . -tao_lmm_rescale_beta - beta factor for rescaling diagonal 2831522df2eSJason Sarich . -tao_lmm_scalar_history - amount of history for scalar scaling 2841522df2eSJason Sarich . -tao_lmm_rescale_history - amount of history for rescaling diagonal 2851522df2eSJason Sarich - -tao_lmm_eps - rejection tolerance 2861522df2eSJason Sarich 2871eb8069cSJason Sarich Level: beginner 2881522df2eSJason Sarich M*/ 289a7e14dcfSSatish Balay #undef __FUNCT__ 290a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_BLMVM" 291728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 292a7e14dcfSSatish Balay { 293a7e14dcfSSatish Balay TAO_BLMVM *blmP; 2948caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 295a7e14dcfSSatish Balay PetscErrorCode ierr; 296a7e14dcfSSatish Balay 297a7e14dcfSSatish Balay PetscFunctionBegin; 298a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BLMVM; 299a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BLMVM; 300a7e14dcfSSatish Balay tao->ops->view = TaoView_BLMVM; 301a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 302a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BLMVM; 303a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_BLMVM; 304a7e14dcfSSatish Balay 3053c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr); 306a9603a14SPatrick Farrell blmP->H0 = NULL; 307a7e14dcfSSatish Balay tao->data = (void*)blmP; 3086552cf8aSJason Sarich 3096552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3106552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 3116552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 312a7e14dcfSSatish Balay 313a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 314a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 315441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 3165d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 317a7e14dcfSSatish Balay PetscFunctionReturn(0); 318a7e14dcfSSatish Balay } 319a7e14dcfSSatish Balay 320a9603a14SPatrick Farrell #undef __FUNCT__ 321a9603a14SPatrick Farrell #define __FUNCT__ "TaoLMVMSetH0" 322a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 323a9603a14SPatrick Farrell { 324a9603a14SPatrick Farrell TAO_LMVM *lmP; 325a9603a14SPatrick Farrell TAO_BLMVM *blmP; 326a9603a14SPatrick Farrell const TaoType type; 327a9603a14SPatrick Farrell PetscBool is_lmvm, is_blmvm; 328a9603a14SPatrick Farrell 329a9603a14SPatrick Farrell PetscErrorCode ierr; 330a9603a14SPatrick Farrell 331a9603a14SPatrick Farrell ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 332a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 333a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 334a9603a14SPatrick Farrell 335a9603a14SPatrick Farrell if (is_lmvm) { 336a9603a14SPatrick Farrell lmP = (TAO_LMVM *)tao->data; 337a9603a14SPatrick Farrell ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 338a9603a14SPatrick Farrell lmP->H0 = H0; 339a9603a14SPatrick Farrell } else if (is_blmvm) { 340a9603a14SPatrick Farrell blmP = (TAO_BLMVM *)tao->data; 341a9603a14SPatrick Farrell ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 342a9603a14SPatrick Farrell blmP->H0 = H0; 34374c66251SPatrick Farrell } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 344a9603a14SPatrick Farrell 345a9603a14SPatrick Farrell PetscFunctionReturn(0); 346a9603a14SPatrick Farrell } 347a9603a14SPatrick Farrell 348a9603a14SPatrick Farrell #undef __FUNCT__ 349a9603a14SPatrick Farrell #define __FUNCT__ "TaoLMVMGetH0" 350a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 351a9603a14SPatrick Farrell { 352a9603a14SPatrick Farrell TAO_LMVM *lmP; 353a9603a14SPatrick Farrell TAO_BLMVM *blmP; 354a9603a14SPatrick Farrell const TaoType type; 355a9603a14SPatrick Farrell PetscBool is_lmvm, is_blmvm; 356a9603a14SPatrick Farrell Mat M; 357a9603a14SPatrick Farrell 358a9603a14SPatrick Farrell PetscErrorCode ierr; 359a9603a14SPatrick Farrell 360a9603a14SPatrick Farrell ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 361a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 362a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 363a9603a14SPatrick Farrell 364a9603a14SPatrick Farrell if (is_lmvm) { 365a9603a14SPatrick Farrell lmP = (TAO_LMVM *)tao->data; 366a9603a14SPatrick Farrell M = lmP->M; 367a9603a14SPatrick Farrell } else if (is_blmvm) { 368a9603a14SPatrick Farrell blmP = (TAO_BLMVM *)tao->data; 369a9603a14SPatrick Farrell M = blmP->M; 37074c66251SPatrick Farrell } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 371a9603a14SPatrick Farrell 372a9603a14SPatrick Farrell ierr = MatLMVMGetH0(M, H0);CHKERRQ(ierr); 373a9603a14SPatrick Farrell PetscFunctionReturn(0); 374a9603a14SPatrick Farrell } 375a9603a14SPatrick Farrell 376a9603a14SPatrick Farrell #undef __FUNCT__ 377a9603a14SPatrick Farrell #define __FUNCT__ "TaoLMVMGetH0KSP" 378a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 379a9603a14SPatrick Farrell { 380a9603a14SPatrick Farrell TAO_LMVM *lmP; 381a9603a14SPatrick Farrell TAO_BLMVM *blmP; 382a9603a14SPatrick Farrell const TaoType type; 383a9603a14SPatrick Farrell PetscBool is_lmvm, is_blmvm; 384a9603a14SPatrick Farrell Mat M; 385a9603a14SPatrick Farrell PetscErrorCode ierr; 386a9603a14SPatrick Farrell 387a9603a14SPatrick Farrell ierr = TaoGetType(tao, &type);CHKERRQ(ierr); 388a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 389a9603a14SPatrick Farrell ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr); 390a9603a14SPatrick Farrell 391a9603a14SPatrick Farrell if (is_lmvm) { 392a9603a14SPatrick Farrell lmP = (TAO_LMVM *)tao->data; 393a9603a14SPatrick Farrell M = lmP->M; 394a9603a14SPatrick Farrell } else if (is_blmvm) { 395a9603a14SPatrick Farrell blmP = (TAO_BLMVM *)tao->data; 396a9603a14SPatrick Farrell M = blmP->M; 39774c66251SPatrick Farrell } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM."); 398a9603a14SPatrick Farrell 399a9603a14SPatrick Farrell ierr = MatLMVMGetH0KSP(M, ksp);CHKERRQ(ierr); 400a9603a14SPatrick Farrell PetscFunctionReturn(0); 401a9603a14SPatrick Farrell } 402