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