1aaa7dc30SBarry Smith #include <../src/tao/pde_constrained/impls/lcl/lcl.h> 2a7e14dcfSSatish Balay static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *); 3a7e14dcfSSatish Balay static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *); 4a7e14dcfSSatish Balay static PetscErrorCode LCLScatter(TAO_LCL *, Vec, Vec, Vec); 5a7e14dcfSSatish Balay static PetscErrorCode LCLGather(TAO_LCL *, Vec, Vec, Vec); 6a7e14dcfSSatish Balay 79371c9d4SSatish Balay static PetscErrorCode TaoDestroy_LCL(Tao tao) { 8a7e14dcfSSatish Balay TAO_LCL *lclP = (TAO_LCL *)tao->data; 9f06e3bfaSBarry Smith 10a7e14dcfSSatish Balay PetscFunctionBegin; 11a7e14dcfSSatish Balay if (tao->setupcalled) { 129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lclP->R)); 139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->lamda)); 149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->lamda0)); 159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->WL)); 169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->W)); 179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->X0)); 189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->G0)); 199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GL)); 209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GAugL)); 219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->dbar)); 229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->U)); 239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->U0)); 249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->V)); 259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->V0)); 269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->V1)); 279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GU)); 289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GV)); 299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GU0)); 309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GV0)); 319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GL_U)); 329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GL_V)); 339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GAugL_U)); 349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GAugL_V)); 359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GL_U0)); 369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GL_V0)); 379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GAugL_U0)); 389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->GAugL_V0)); 399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->DU)); 409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->DV)); 419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->WU)); 429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->WV)); 439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->g1)); 449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->g2)); 459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->con1)); 46a7e14dcfSSatish Balay 479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->r)); 489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lclP->s)); 49a7e14dcfSSatish Balay 509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tao->state_is)); 519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tao->design_is)); 52a7e14dcfSSatish Balay 539566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&lclP->state_scatter)); 549566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&lclP->design_scatter)); 55a7e14dcfSSatish Balay } 569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lclP->R)); 57a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 589566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 59a7e14dcfSSatish Balay PetscFunctionReturn(0); 60a7e14dcfSSatish Balay } 61a7e14dcfSSatish Balay 629371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_LCL(Tao tao, PetscOptionItems *PetscOptionsObject) { 63a7e14dcfSSatish Balay TAO_LCL *lclP = (TAO_LCL *)tao->data; 64a7e14dcfSSatish Balay 65a7e14dcfSSatish Balay PetscFunctionBegin; 66d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Linearly-Constrained Augmented Lagrangian Method for PDE-constrained optimization"); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_eps1", "epsilon 1 tolerance", "", lclP->eps1, &lclP->eps1, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_eps2", "epsilon 2 tolerance", "", lclP->eps2, &lclP->eps2, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_rho0", "init value for rho", "", lclP->rho0, &lclP->rho0, NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_rhomax", "max value for rho", "", lclP->rhomax, &lclP->rhomax, NULL)); 71a7e14dcfSSatish Balay lclP->phase2_niter = 1; 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_lcl_phase2_niter", "Number of phase 2 iterations in LCL algorithm", "", lclP->phase2_niter, &lclP->phase2_niter, NULL)); 73a7e14dcfSSatish Balay lclP->verbose = PETSC_FALSE; 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_lcl_verbose", "Print verbose output", "", lclP->verbose, &lclP->verbose, NULL)); 75a7e14dcfSSatish Balay lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4; 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_tola", "Tolerance for first forward solve", "", lclP->tau[0], &lclP->tau[0], NULL)); 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_tolb", "Tolerance for first adjoint solve", "", lclP->tau[1], &lclP->tau[1], NULL)); 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_tolc", "Tolerance for second forward solve", "", lclP->tau[2], &lclP->tau[2], NULL)); 799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_lcl_told", "Tolerance for second adjoint solve", "", lclP->tau[3], &lclP->tau[3], NULL)); 80d0609cedSBarry Smith PetscOptionsHeadEnd(); 819566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 829566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(lclP->R)); 83a7e14dcfSSatish Balay PetscFunctionReturn(0); 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay 869371c9d4SSatish Balay static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer) { 87a7e14dcfSSatish Balay return 0; 88a7e14dcfSSatish Balay } 89a7e14dcfSSatish Balay 909371c9d4SSatish Balay static PetscErrorCode TaoSetup_LCL(Tao tao) { 91a7e14dcfSSatish Balay TAO_LCL *lclP = (TAO_LCL *)tao->data; 92a7e14dcfSSatish Balay PetscInt lo, hi, nlocalstate, nlocaldesign; 93a7e14dcfSSatish Balay IS is_state, is_design; 94f06e3bfaSBarry Smith 95a7e14dcfSSatish Balay PetscFunctionBegin; 963c859ba3SBarry Smith PetscCheck(tao->state_is, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "LCL Solver requires an initial state index set -- use TaoSetStateIS()"); 979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &lclP->W)); 1009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &lclP->X0)); 1019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &lclP->G0)); 1029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &lclP->GL)); 1039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &lclP->GAugL)); 104a7e14dcfSSatish Balay 1059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints, &lclP->lamda)); 1069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints, &lclP->WL)); 1079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints, &lclP->lamda0)); 1089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints, &lclP->con1)); 109a7e14dcfSSatish Balay 1109566063dSJacob Faibussowitsch PetscCall(VecSet(lclP->lamda, 0.0)); 111a7e14dcfSSatish Balay 1129566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &lclP->n)); 1139566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints, &lclP->m)); 114a7e14dcfSSatish Balay 1159566063dSJacob Faibussowitsch PetscCall(VecCreate(((PetscObject)tao)->comm, &lclP->U)); 1169566063dSJacob Faibussowitsch PetscCall(VecCreate(((PetscObject)tao)->comm, &lclP->V)); 1179566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(tao->state_is, &nlocalstate)); 1189566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(tao->design_is, &nlocaldesign)); 1199566063dSJacob Faibussowitsch PetscCall(VecSetSizes(lclP->U, nlocalstate, lclP->m)); 1209566063dSJacob Faibussowitsch PetscCall(VecSetSizes(lclP->V, nlocaldesign, lclP->n - lclP->m)); 1219566063dSJacob Faibussowitsch PetscCall(VecSetType(lclP->U, ((PetscObject)(tao->solution))->type_name)); 1229566063dSJacob Faibussowitsch PetscCall(VecSetType(lclP->V, ((PetscObject)(tao->solution))->type_name)); 1239566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(lclP->U)); 1249566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(lclP->V)); 1259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->DU)); 1269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->U0)); 1279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->GU)); 1289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->GU0)); 1299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->GAugL_U)); 1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->GL_U)); 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->GAugL_U0)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->GL_U0)); 1339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->WU)); 1349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->U, &lclP->r)); 1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->V0)); 1369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->V1)); 1379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->DV)); 1389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->s)); 1399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->GV)); 1409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->GV0)); 1419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->dbar)); 1429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->GAugL_V)); 1439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->GL_V)); 1449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->GAugL_V0)); 1459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->GL_V0)); 1469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->WV)); 1479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->g1)); 1489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lclP->V, &lclP->g2)); 149a7e14dcfSSatish Balay 150a7e14dcfSSatish Balay /* create scatters for state, design subvecs */ 1519566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(lclP->U, &lo, &hi)); 1529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(((PetscObject)lclP->U)->comm, hi - lo, lo, 1, &is_state)); 1539566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(lclP->V, &lo, &hi)); 154a7e14dcfSSatish Balay if (0) { 155a7e14dcfSSatish Balay PetscInt sizeU, sizeV; 1569566063dSJacob Faibussowitsch PetscCall(VecGetSize(lclP->U, &sizeU)); 1579566063dSJacob Faibussowitsch PetscCall(VecGetSize(lclP->V, &sizeV)); 15863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "size(U)=%" PetscInt_FMT ", size(V)=%" PetscInt_FMT "\n", sizeU, sizeV)); 159a7e14dcfSSatish Balay } 1609566063dSJacob Faibussowitsch PetscCall(ISCreateStride(((PetscObject)lclP->V)->comm, hi - lo, lo, 1, &is_design)); 1619566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->solution, tao->state_is, lclP->U, is_state, &lclP->state_scatter)); 1629566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->solution, tao->design_is, lclP->V, is_design, &lclP->design_scatter)); 1639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_state)); 1649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_design)); 165a7e14dcfSSatish Balay PetscFunctionReturn(0); 166a7e14dcfSSatish Balay } 167a7e14dcfSSatish Balay 1689371c9d4SSatish Balay static PetscErrorCode TaoSolve_LCL(Tao tao) { 169a7e14dcfSSatish Balay TAO_LCL *lclP = (TAO_LCL *)tao->data; 1708931d482SJason Sarich PetscInt phase2_iter, nlocal, its; 171e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 172a7e14dcfSSatish Balay PetscReal step = 1.0, f, descent, aldescent; 173a7e14dcfSSatish Balay PetscReal cnorm, mnorm; 174a7e14dcfSSatish Balay PetscReal adec, r2, rGL_U, rWU; 175a7e14dcfSSatish Balay PetscBool set, pset, flag, pflag, symmetric; 176a7e14dcfSSatish Balay 177f06e3bfaSBarry Smith PetscFunctionBegin; 178a7e14dcfSSatish Balay lclP->rho = lclP->rho0; 1799566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(lclP->U, &nlocal)); 1809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(lclP->V, &nlocal)); 1819566063dSJacob Faibussowitsch PetscCall(MatSetSizes(lclP->R, nlocal, nlocal, lclP->n - lclP->m, lclP->n - lclP->m)); 1829566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(lclP->R, lclP->V, lclP->V)); 183a7e14dcfSSatish Balay lclP->recompute_jacobian_flag = PETSC_TRUE; 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay /* Scatter to U,V */ 1869566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V)); 187a7e14dcfSSatish Balay 188a7e14dcfSSatish Balay /* Evaluate Function, Gradient, Constraints, and Jacobian */ 1899566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 1909566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianState(tao, tao->solution, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv)); 1919566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianDesign(tao, tao->solution, tao->jacobian_design)); 1929566063dSJacob Faibussowitsch PetscCall(TaoComputeConstraints(tao, tao->solution, tao->constraints)); 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay /* Scatter gradient to GU,GV */ 1959566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV)); 196a7e14dcfSSatish Balay 197a7e14dcfSSatish Balay /* Evaluate Lagrangian function and gradient */ 198a7e14dcfSSatish Balay /* p0 */ 1999566063dSJacob Faibussowitsch PetscCall(VecSet(lclP->lamda, 0.0)); /* Initial guess in CG */ 2009566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag)); 201a7e14dcfSSatish Balay if (tao->jacobian_state_pre) { 2029566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag)); 203a7e14dcfSSatish Balay } else { 204a7e14dcfSSatish Balay pset = pflag = PETSC_TRUE; 205a7e14dcfSSatish Balay } 206f06e3bfaSBarry Smith if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 207f06e3bfaSBarry Smith else symmetric = PETSC_FALSE; 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay lclP->solve_type = LCL_ADJOINT2; 210a7e14dcfSSatish Balay if (tao->jacobian_state_inv) { 211a7e14dcfSSatish Balay if (symmetric) { 2129371c9d4SSatish Balay PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); 2139371c9d4SSatish Balay } else { 2149566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay } else { 2179566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre)); 218a7e14dcfSSatish Balay if (symmetric) { 2199566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, lclP->GU, lclP->lamda)); 220a7e14dcfSSatish Balay } else { 2219566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda)); 222a7e14dcfSSatish Balay } 2239566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 224a7e14dcfSSatish Balay tao->ksp_its += its; 225ae93cb3cSJason Sarich tao->ksp_tot_its += its; 226a7e14dcfSSatish Balay } 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->lamda, lclP->lamda0)); 2289566063dSJacob Faibussowitsch PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao)); 229a7e14dcfSSatish Balay 2309566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GL, lclP->GL_U, lclP->GL_V)); 2319566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V)); 232a7e14dcfSSatish Balay 233a7e14dcfSSatish Balay /* Evaluate constraint norm */ 2349566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm)); 2359566063dSJacob Faibussowitsch PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm)); 236a7e14dcfSSatish Balay 237a7e14dcfSSatish Balay /* Monitor convergence */ 238763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 2399566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, mnorm, cnorm, tao->ksp_its)); 2409566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, mnorm, cnorm, step)); 241dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 242a7e14dcfSSatish Balay 243763847b4SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 244e1e80dc8SAlp Dener /* Call general purpose update function */ 245dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 246ae93cb3cSJason Sarich tao->ksp_its = 0; 247a7e14dcfSSatish Balay /* Compute a descent direction for the linearly constrained subproblem 248a7e14dcfSSatish Balay minimize f(u+du, v+dv) 249a7e14dcfSSatish Balay s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */ 250a7e14dcfSSatish Balay 251a7e14dcfSSatish Balay /* Store the points around the linearization */ 2529566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->U, lclP->U0)); 2539566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->V, lclP->V0)); 2549566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GU, lclP->GU0)); 2559566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GV, lclP->GV0)); 2569566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GAugL_U, lclP->GAugL_U0)); 2579566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GAugL_V, lclP->GAugL_V0)); 2589566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GL_U, lclP->GL_U0)); 2599566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GL_V, lclP->GL_V0)); 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay lclP->aug0 = lclP->aug; 262a7e14dcfSSatish Balay lclP->lgn0 = lclP->lgn; 263a7e14dcfSSatish Balay 264a7e14dcfSSatish Balay /* Given the design variables, we need to project the current iterate 265a7e14dcfSSatish Balay onto the linearized constraint. We choose to fix the design variables 266a7e14dcfSSatish Balay and solve the linear system for the state variables. The resulting 267a7e14dcfSSatish Balay point is the Newton direction */ 268a7e14dcfSSatish Balay 269a7e14dcfSSatish Balay /* Solve r = A\con */ 270a7e14dcfSSatish Balay lclP->solve_type = LCL_FORWARD1; 2719566063dSJacob Faibussowitsch PetscCall(VecSet(lclP->r, 0.0)); /* Initial guess in CG */ 272a7e14dcfSSatish Balay 273a7e14dcfSSatish Balay if (tao->jacobian_state_inv) { 2749566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r)); 275a7e14dcfSSatish Balay } else { 2769566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre)); 2779566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tao->constraints, lclP->r)); 2789566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 279a7e14dcfSSatish Balay tao->ksp_its += its; 280ae93cb3cSJason Sarich tao->ksp_tot_its += tao->ksp_its; 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay 283a7e14dcfSSatish Balay /* Set design step direction dv to zero */ 2849566063dSJacob Faibussowitsch PetscCall(VecSet(lclP->s, 0.0)); 285a7e14dcfSSatish Balay 286a7e14dcfSSatish Balay /* 287a7e14dcfSSatish Balay Check sufficient descent for constraint merit function .5*||con||^2 288a7e14dcfSSatish Balay con' Ak r >= eps1 ||r||^(2+eps2) 289a7e14dcfSSatish Balay */ 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay /* Compute WU= Ak' * con */ 292a7e14dcfSSatish Balay if (symmetric) { 2939566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state, tao->constraints, lclP->WU)); 294a7e14dcfSSatish Balay } else { 2959566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_state, tao->constraints, lclP->WU)); 296a7e14dcfSSatish Balay } 297a7e14dcfSSatish Balay /* Compute r * Ak' * con */ 2989566063dSJacob Faibussowitsch PetscCall(VecDot(lclP->r, lclP->WU, &rWU)); 299a7e14dcfSSatish Balay 300a7e14dcfSSatish Balay /* compute ||r||^(2+eps2) */ 3019566063dSJacob Faibussowitsch PetscCall(VecNorm(lclP->r, NORM_2, &r2)); 302a7e14dcfSSatish Balay r2 = PetscPowScalar(r2, 2.0 + lclP->eps2); 303a7e14dcfSSatish Balay adec = lclP->eps1 * r2; 304a7e14dcfSSatish Balay 305a7e14dcfSSatish Balay if (rWU < adec) { 3069566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Newton direction not descent for constraint, feasibility phase required\n")); 30748a46eb9SPierre Jolivet if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Newton direction not descent for constraint: %g -- using steepest descent\n", (double)descent)); 308a7e14dcfSSatish Balay 3099566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Using steepest descent direction instead.\n")); 3109566063dSJacob Faibussowitsch PetscCall(VecSet(lclP->r, 0.0)); 3119566063dSJacob Faibussowitsch PetscCall(VecAXPY(lclP->r, -1.0, lclP->WU)); 3129566063dSJacob Faibussowitsch PetscCall(VecDot(lclP->r, lclP->r, &rWU)); 3139566063dSJacob Faibussowitsch PetscCall(VecNorm(lclP->r, NORM_2, &r2)); 314a7e14dcfSSatish Balay r2 = PetscPowScalar(r2, 2.0 + lclP->eps2); 3159566063dSJacob Faibussowitsch PetscCall(VecDot(lclP->r, lclP->GAugL_U, &descent)); 316a7e14dcfSSatish Balay adec = lclP->eps1 * r2; 317a7e14dcfSSatish Balay } 318a7e14dcfSSatish Balay 319a7e14dcfSSatish Balay /* 320a7e14dcfSSatish Balay Check descent for aug. lagrangian 321a7e14dcfSSatish Balay r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2) 322a7e14dcfSSatish Balay GL_U = GUk - Ak'*yk 323a7e14dcfSSatish Balay WU = Ak'*con 324a7e14dcfSSatish Balay adec=eps1||r||^(2+eps2) 325a7e14dcfSSatish Balay 326a7e14dcfSSatish Balay ==> 327a7e14dcfSSatish Balay Check r'GL_U - rho*r'WU <= adec 328a7e14dcfSSatish Balay */ 329a7e14dcfSSatish Balay 3309566063dSJacob Faibussowitsch PetscCall(VecDot(lclP->r, lclP->GL_U, &rGL_U)); 331a7e14dcfSSatish Balay aldescent = rGL_U - lclP->rho * rWU; 332a7e14dcfSSatish Balay if (aldescent > -adec) { 33348a46eb9SPierre Jolivet if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Newton direction not descent for augmented Lagrangian: %g", (double)aldescent)); 3349566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Newton direction not descent for augmented Lagrangian: %g\n", (double)aldescent)); 335a7e14dcfSSatish Balay lclP->rho = (rGL_U - adec) / rWU; 336a7e14dcfSSatish Balay if (lclP->rho > lclP->rhomax) { 337a7e14dcfSSatish Balay lclP->rho = lclP->rhomax; 33898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "rho=%g > rhomax, case not implemented. Increase rhomax (-tao_lcl_rhomax)", (double)lclP->rho); 339a7e14dcfSSatish Balay } 34048a46eb9SPierre Jolivet if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Increasing penalty parameter to %g\n", (double)lclP->rho)); 3419566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, " Increasing penalty parameter to %g\n", (double)lclP->rho)); 342a7e14dcfSSatish Balay } 343a7e14dcfSSatish Balay 3449566063dSJacob Faibussowitsch PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao)); 3459566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GL, lclP->GL_U, lclP->GL_V)); 3469566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V)); 347a7e14dcfSSatish Balay 348a7e14dcfSSatish Balay /* We now minimize the augmented Lagrangian along the Newton direction */ 3499566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->r, -1.0)); 3509566063dSJacob Faibussowitsch PetscCall(LCLGather(lclP, lclP->r, lclP->s, tao->stepdirection)); 3519566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->r, -1.0)); 3529566063dSJacob Faibussowitsch PetscCall(LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL)); 3539566063dSJacob Faibussowitsch PetscCall(LCLGather(lclP, lclP->U0, lclP->V0, lclP->X0)); 354a7e14dcfSSatish Balay 355a7e14dcfSSatish Balay lclP->recompute_jacobian_flag = PETSC_TRUE; 356a7e14dcfSSatish Balay 3579566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 3589566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT)); 3599566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 3609566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason)); 36148a46eb9SPierre Jolivet if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steplength = %g\n", (double)step)); 362a7e14dcfSSatish Balay 3639566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V)); 3649566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 3659566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV)); 366a7e14dcfSSatish Balay 3679566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V)); 368a7e14dcfSSatish Balay 369a7e14dcfSSatish Balay /* Check convergence */ 3709566063dSJacob Faibussowitsch PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm)); 3719566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm)); 3729566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, mnorm, cnorm, tao->ksp_its)); 3739566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, mnorm, cnorm, step)); 374dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 375ad540459SPierre Jolivet if (tao->reason != TAO_CONTINUE_ITERATING) break; 376a7e14dcfSSatish Balay 377a7e14dcfSSatish Balay /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */ 378a7e14dcfSSatish Balay for (phase2_iter = 0; phase2_iter < lclP->phase2_niter; phase2_iter++) { 379a7e14dcfSSatish Balay /* We now minimize the objective function starting from the fraction of 380a7e14dcfSSatish Balay the Newton point accepted by applying one step of a reduced-space 381a7e14dcfSSatish Balay method. The optimization problem is: 382a7e14dcfSSatish Balay 383a7e14dcfSSatish Balay minimize f(u+du, v+dv) 384a7e14dcfSSatish Balay s. t. A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0) 385a7e14dcfSSatish Balay 386a7e14dcfSSatish Balay In particular, we have that 387a7e14dcfSSatish Balay du = -inv(A)*(Bdv + alpha g) */ 388a7e14dcfSSatish Balay 3899566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianState(tao, lclP->X0, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv)); 3909566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianDesign(tao, lclP->X0, tao->jacobian_design)); 391a7e14dcfSSatish Balay 392a7e14dcfSSatish Balay /* Store V and constraints */ 3939566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->V, lclP->V1)); 3949566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->constraints, lclP->con1)); 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay /* Compute multipliers */ 397a7e14dcfSSatish Balay /* p1 */ 3989566063dSJacob Faibussowitsch PetscCall(VecSet(lclP->lamda, 0.0)); /* Initial guess in CG */ 399a7e14dcfSSatish Balay lclP->solve_type = LCL_ADJOINT1; 4009566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag)); 401a7e14dcfSSatish Balay if (tao->jacobian_state_pre) { 4029566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag)); 403a7e14dcfSSatish Balay } else { 404a7e14dcfSSatish Balay pset = pflag = PETSC_TRUE; 405a7e14dcfSSatish Balay } 406f06e3bfaSBarry Smith if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 407f06e3bfaSBarry Smith else symmetric = PETSC_FALSE; 408a7e14dcfSSatish Balay 409a7e14dcfSSatish Balay if (tao->jacobian_state_inv) { 410a7e14dcfSSatish Balay if (symmetric) { 4119566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda)); 412a7e14dcfSSatish Balay } else { 4139566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda)); 414a7e14dcfSSatish Balay } 415a7e14dcfSSatish Balay } else { 416a7e14dcfSSatish Balay if (symmetric) { 4179566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, lclP->GAugL_U, lclP->lamda)); 418a7e14dcfSSatish Balay } else { 4199566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(tao->ksp, lclP->GAugL_U, lclP->lamda)); 420a7e14dcfSSatish Balay } 4219566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 422a7e14dcfSSatish Balay tao->ksp_its += its; 423ae93cb3cSJason Sarich tao->ksp_tot_its += its; 424a7e14dcfSSatish Balay } 4259566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_design, lclP->lamda, lclP->g1)); 4269566063dSJacob Faibussowitsch PetscCall(VecAXPY(lclP->g1, -1.0, lclP->GAugL_V)); 427a7e14dcfSSatish Balay 428a7e14dcfSSatish Balay /* Compute the limited-memory quasi-newton direction */ 4298931d482SJason Sarich if (tao->niter > 0) { 4309566063dSJacob Faibussowitsch PetscCall(MatSolve(lclP->R, lclP->g1, lclP->s)); 4319566063dSJacob Faibussowitsch PetscCall(VecDot(lclP->s, lclP->g1, &descent)); 4320583ad50SJason Sarich if (descent <= 0) { 43348a46eb9SPierre Jolivet if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Reduced-space direction not descent: %g\n", (double)descent)); 4349566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->g1, lclP->s)); 435a7e14dcfSSatish Balay } 436a7e14dcfSSatish Balay } else { 4379566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->g1, lclP->s)); 438a7e14dcfSSatish Balay } 4399566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->g1, -1.0)); 440a7e14dcfSSatish Balay 441a7e14dcfSSatish Balay /* Recover the full space direction */ 4429566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_design, lclP->s, lclP->WU)); 4439566063dSJacob Faibussowitsch /* PetscCall(VecSet(lclP->r,0.0)); */ /* Initial guess in CG */ 444a7e14dcfSSatish Balay lclP->solve_type = LCL_FORWARD2; 445a7e14dcfSSatish Balay if (tao->jacobian_state_inv) { 4469566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state_inv, lclP->WU, lclP->r)); 447a7e14dcfSSatish Balay } else { 4489566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, lclP->WU, lclP->r)); 4499566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 450a7e14dcfSSatish Balay tao->ksp_its += its; 451ae93cb3cSJason Sarich tao->ksp_tot_its += its; 452a7e14dcfSSatish Balay } 453a7e14dcfSSatish Balay 454a7e14dcfSSatish Balay /* We now minimize the augmented Lagrangian along the direction -r,s */ 4559566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->r, -1.0)); 4569566063dSJacob Faibussowitsch PetscCall(LCLGather(lclP, lclP->r, lclP->s, tao->stepdirection)); 4579566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->r, -1.0)); 458a7e14dcfSSatish Balay lclP->recompute_jacobian_flag = PETSC_TRUE; 459a7e14dcfSSatish Balay 4609566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 4619566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT)); 4629566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 4639566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason)); 46448a46eb9SPierre Jolivet if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Reduced-space steplength = %g\n", (double)step)); 465a7e14dcfSSatish Balay 4669566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V)); 4679566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GL, lclP->GL_U, lclP->GL_V)); 4689566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V)); 4699566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 4709566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV)); 471a7e14dcfSSatish Balay 472a7e14dcfSSatish Balay /* Compute the reduced gradient at the new point */ 473a7e14dcfSSatish Balay 4749566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianState(tao, lclP->X0, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv)); 4759566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianDesign(tao, lclP->X0, tao->jacobian_design)); 476a7e14dcfSSatish Balay 477a7e14dcfSSatish Balay /* p2 */ 478a7e14dcfSSatish Balay /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */ 479a7e14dcfSSatish Balay if (phase2_iter == 0) { 4809566063dSJacob Faibussowitsch PetscCall(VecWAXPY(lclP->lamda, -lclP->rho, lclP->con1, lclP->lamda0)); 481f06e3bfaSBarry Smith } else { 4829566063dSJacob Faibussowitsch PetscCall(VecAXPY(lclP->lamda, -lclP->rho, tao->constraints)); 483a7e14dcfSSatish Balay } 484a7e14dcfSSatish Balay 4859566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag)); 486a7e14dcfSSatish Balay if (tao->jacobian_state_pre) { 4879566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag)); 488a7e14dcfSSatish Balay } else { 489a7e14dcfSSatish Balay pset = pflag = PETSC_TRUE; 490a7e14dcfSSatish Balay } 491f06e3bfaSBarry Smith if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 492f06e3bfaSBarry Smith else symmetric = PETSC_FALSE; 493a7e14dcfSSatish Balay 494a7e14dcfSSatish Balay lclP->solve_type = LCL_ADJOINT2; 495a7e14dcfSSatish Balay if (tao->jacobian_state_inv) { 496a7e14dcfSSatish Balay if (symmetric) { 4979566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); 498a7e14dcfSSatish Balay } else { 4999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); 500a7e14dcfSSatish Balay } 501a7e14dcfSSatish Balay } else { 502a7e14dcfSSatish Balay if (symmetric) { 5039566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, lclP->GU, lclP->lamda)); 504a7e14dcfSSatish Balay } else { 5059566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda)); 506a7e14dcfSSatish Balay } 5079566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 508a7e14dcfSSatish Balay tao->ksp_its += its; 5092d9aa51bSJason Sarich tao->ksp_tot_its += its; 510a7e14dcfSSatish Balay } 511a7e14dcfSSatish Balay 5129566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_design, lclP->lamda, lclP->g2)); 5139566063dSJacob Faibussowitsch PetscCall(VecAXPY(lclP->g2, -1.0, lclP->GV)); 514a7e14dcfSSatish Balay 5159566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->g2, -1.0)); 516a7e14dcfSSatish Balay 517a7e14dcfSSatish Balay /* Update the quasi-newton approximation */ 5189566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(lclP->R, lclP->V, lclP->g2)); 519*750b007cSBarry Smith /* Use "-tao_ls_type gpcg -tao_ls_ftol 0 -tao_lmm_broyden_phi 0.0 -tao_lmm_scale_type scalar" to obtain agreement with MATLAB code */ 520a7e14dcfSSatish Balay } 521a7e14dcfSSatish Balay 5229566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->lamda, lclP->lamda0)); 523a7e14dcfSSatish Balay 524a7e14dcfSSatish Balay /* Evaluate Function, Gradient, Constraints, and Jacobian */ 5259566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 5269566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V)); 5279566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV)); 528a7e14dcfSSatish Balay 5299566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianState(tao, tao->solution, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv)); 5309566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianDesign(tao, tao->solution, tao->jacobian_design)); 5319566063dSJacob Faibussowitsch PetscCall(TaoComputeConstraints(tao, tao->solution, tao->constraints)); 532a7e14dcfSSatish Balay 5339566063dSJacob Faibussowitsch PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao)); 534a7e14dcfSSatish Balay 5359566063dSJacob Faibussowitsch PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm)); 536a7e14dcfSSatish Balay 537a7e14dcfSSatish Balay /* Evaluate constraint norm */ 5389566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm)); 539a7e14dcfSSatish Balay 540a7e14dcfSSatish Balay /* Monitor convergence */ 5418931d482SJason Sarich tao->niter++; 5429566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, mnorm, cnorm, tao->ksp_its)); 5439566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, mnorm, cnorm, step)); 544dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 545a7e14dcfSSatish Balay } 546a7e14dcfSSatish Balay PetscFunctionReturn(0); 547a7e14dcfSSatish Balay } 548a7e14dcfSSatish Balay 5491522df2eSJason Sarich /*MC 5501522df2eSJason Sarich TAOLCL - linearly constrained lagrangian method for pde-constrained optimization 5511522df2eSJason Sarich 5521522df2eSJason Sarich + -tao_lcl_eps1 - epsilon 1 tolerance 553d0609cedSBarry Smith . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL); 554d0609cedSBarry Smith . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL); 555d0609cedSBarry Smith . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL); 5561522df2eSJason Sarich . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm 5571522df2eSJason Sarich . -tao_lcl_verbose - Print verbose output if True 5581522df2eSJason Sarich . -tao_lcl_tola - Tolerance for first forward solve 5591522df2eSJason Sarich . -tao_lcl_tolb - Tolerance for first adjoint solve 5601522df2eSJason Sarich . -tao_lcl_tolc - Tolerance for second forward solve 5611522df2eSJason Sarich - -tao_lcl_told - Tolerance for second adjoint solve 5621522df2eSJason Sarich 5631eb8069cSJason Sarich Level: beginner 5641522df2eSJason Sarich M*/ 5659371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao) { 566a7e14dcfSSatish Balay TAO_LCL *lclP; 5678caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 568a7e14dcfSSatish Balay 569f06e3bfaSBarry Smith PetscFunctionBegin; 570a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_LCL; 571a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LCL; 572a7e14dcfSSatish Balay tao->ops->view = TaoView_LCL; 573a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LCL; 574a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LCL; 5754dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lclP)); 576a7e14dcfSSatish Balay tao->data = (void *)lclP; 577a7e14dcfSSatish Balay 5786552cf8aSJason Sarich /* Override default settings (unless already changed) */ 5796552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 200; 5806552cf8aSJason Sarich if (!tao->catol_changed) tao->catol = 1.0e-4; 5816552cf8aSJason Sarich if (!tao->gatol_changed) tao->gttol = 1.0e-4; 5826552cf8aSJason Sarich if (!tao->grtol_changed) tao->gttol = 1.0e-4; 5836552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 1.0e-4; 584a7e14dcfSSatish Balay lclP->rho0 = 1.0e-4; 585a7e14dcfSSatish Balay lclP->rhomax = 1e5; 586a7e14dcfSSatish Balay lclP->eps1 = 1.0e-8; 587a7e14dcfSSatish Balay lclP->eps2 = 0.0; 588a7e14dcfSSatish Balay lclP->solve_type = 2; 589a7e14dcfSSatish Balay lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4; 5909566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 5919566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 5929566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 5939566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix)); 594a7e14dcfSSatish Balay 5959566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, LCLComputeAugmentedLagrangianAndGradient, tao)); 5969566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 5979566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 5989566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 5999566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 6000c51296cSAlp Dener 6019566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &lclP->R)); 6029566063dSJacob Faibussowitsch PetscCall(MatSetType(lclP->R, MATLMVMBFGS)); 603a7e14dcfSSatish Balay PetscFunctionReturn(0); 604a7e14dcfSSatish Balay } 605a7e14dcfSSatish Balay 6069371c9d4SSatish Balay static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) { 607441846f8SBarry Smith Tao tao = (Tao)ptr; 608a7e14dcfSSatish Balay TAO_LCL *lclP = (TAO_LCL *)tao->data; 609a7e14dcfSSatish Balay PetscBool set, pset, flag, pflag, symmetric; 610a7e14dcfSSatish Balay PetscReal cdotl; 611a7e14dcfSSatish Balay 612a7e14dcfSSatish Balay PetscFunctionBegin; 6139566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, X, f, G)); 6149566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, G, lclP->GU, lclP->GV)); 615a7e14dcfSSatish Balay if (lclP->recompute_jacobian_flag) { 6169566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianState(tao, X, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv)); 6179566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianDesign(tao, X, tao->jacobian_design)); 618a7e14dcfSSatish Balay } 6199566063dSJacob Faibussowitsch PetscCall(TaoComputeConstraints(tao, X, tao->constraints)); 6209566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag)); 621a7e14dcfSSatish Balay if (tao->jacobian_state_pre) { 6229566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag)); 623a7e14dcfSSatish Balay } else { 624a7e14dcfSSatish Balay pset = pflag = PETSC_TRUE; 625a7e14dcfSSatish Balay } 626f06e3bfaSBarry Smith if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 627f06e3bfaSBarry Smith else symmetric = PETSC_FALSE; 628a7e14dcfSSatish Balay 6299566063dSJacob Faibussowitsch PetscCall(VecDot(lclP->lamda0, tao->constraints, &cdotl)); 630a7e14dcfSSatish Balay lclP->lgn = *f - cdotl; 631a7e14dcfSSatish Balay 632a7e14dcfSSatish Balay /* Gradient of Lagrangian GL = G - J' * lamda */ 633a7e14dcfSSatish Balay /* WU = A' * WL 634a7e14dcfSSatish Balay WV = B' * WL */ 635a7e14dcfSSatish Balay if (symmetric) { 6369566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state, lclP->lamda0, lclP->GL_U)); 637a7e14dcfSSatish Balay } else { 6389566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_state, lclP->lamda0, lclP->GL_U)); 639a7e14dcfSSatish Balay } 6409566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_design, lclP->lamda0, lclP->GL_V)); 6419566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->GL_U, -1.0)); 6429566063dSJacob Faibussowitsch PetscCall(VecScale(lclP->GL_V, -1.0)); 6439566063dSJacob Faibussowitsch PetscCall(VecAXPY(lclP->GL_U, 1.0, lclP->GU)); 6449566063dSJacob Faibussowitsch PetscCall(VecAXPY(lclP->GL_V, 1.0, lclP->GV)); 6459566063dSJacob Faibussowitsch PetscCall(LCLGather(lclP, lclP->GL_U, lclP->GL_V, lclP->GL)); 646a7e14dcfSSatish Balay 647a7e14dcfSSatish Balay f[0] = lclP->lgn; 6489566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GL, G)); 649a7e14dcfSSatish Balay PetscFunctionReturn(0); 650a7e14dcfSSatish Balay } 651a7e14dcfSSatish Balay 6529371c9d4SSatish Balay static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) { 653441846f8SBarry Smith Tao tao = (Tao)ptr; 654a7e14dcfSSatish Balay TAO_LCL *lclP = (TAO_LCL *)tao->data; 655a7e14dcfSSatish Balay PetscReal con2; 656a7e14dcfSSatish Balay PetscBool flag, pflag, set, pset, symmetric; 657a7e14dcfSSatish Balay 658a7e14dcfSSatish Balay PetscFunctionBegin; 6599566063dSJacob Faibussowitsch PetscCall(LCLComputeLagrangianAndGradient(tao->linesearch, X, f, G, tao)); 6609566063dSJacob Faibussowitsch PetscCall(LCLScatter(lclP, G, lclP->GL_U, lclP->GL_V)); 6619566063dSJacob Faibussowitsch PetscCall(VecDot(tao->constraints, tao->constraints, &con2)); 662a7e14dcfSSatish Balay lclP->aug = lclP->lgn + 0.5 * lclP->rho * con2; 663a7e14dcfSSatish Balay 664a7e14dcfSSatish Balay /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */ 665a7e14dcfSSatish Balay /* WU = A' * c 666a7e14dcfSSatish Balay WV = B' * c */ 6679566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag)); 668a7e14dcfSSatish Balay if (tao->jacobian_state_pre) { 6699566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag)); 670a7e14dcfSSatish Balay } else { 671a7e14dcfSSatish Balay pset = pflag = PETSC_TRUE; 672a7e14dcfSSatish Balay } 673f06e3bfaSBarry Smith if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 674f06e3bfaSBarry Smith else symmetric = PETSC_FALSE; 675a7e14dcfSSatish Balay 676a7e14dcfSSatish Balay if (symmetric) { 6779566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian_state, tao->constraints, lclP->GAugL_U)); 678a7e14dcfSSatish Balay } else { 6799566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_state, tao->constraints, lclP->GAugL_U)); 680a7e14dcfSSatish Balay } 681a7e14dcfSSatish Balay 6829566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_design, tao->constraints, lclP->GAugL_V)); 6839566063dSJacob Faibussowitsch PetscCall(VecAYPX(lclP->GAugL_U, lclP->rho, lclP->GL_U)); 6849566063dSJacob Faibussowitsch PetscCall(VecAYPX(lclP->GAugL_V, lclP->rho, lclP->GL_V)); 6859566063dSJacob Faibussowitsch PetscCall(LCLGather(lclP, lclP->GAugL_U, lclP->GAugL_V, lclP->GAugL)); 686a7e14dcfSSatish Balay 687a7e14dcfSSatish Balay f[0] = lclP->aug; 6889566063dSJacob Faibussowitsch PetscCall(VecCopy(lclP->GAugL, G)); 689a7e14dcfSSatish Balay PetscFunctionReturn(0); 690a7e14dcfSSatish Balay } 691a7e14dcfSSatish Balay 6929371c9d4SSatish Balay PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x) { 693a7e14dcfSSatish Balay PetscFunctionBegin; 6949566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE)); 6959566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE)); 6969566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE)); 6979566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE)); 698a7e14dcfSSatish Balay PetscFunctionReturn(0); 699a7e14dcfSSatish Balay } 7009371c9d4SSatish Balay PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v) { 701a7e14dcfSSatish Balay PetscFunctionBegin; 7029566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD)); 7039566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD)); 7049566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD)); 7059566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD)); 706a7e14dcfSSatish Balay PetscFunctionReturn(0); 707a7e14dcfSSatish Balay } 708