xref: /petsc/src/tao/pde_constrained/impls/lcl/lcl.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
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 
TaoDestroy_LCL(Tao tao)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_LCL(Tao tao)
8d71ae5a4SJacob Faibussowitsch {
9a7e14dcfSSatish Balay   TAO_LCL *lclP = (TAO_LCL *)tao->data;
10f06e3bfaSBarry Smith 
11a7e14dcfSSatish Balay   PetscFunctionBegin;
12a7e14dcfSSatish Balay   if (tao->setupcalled) {
139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&lclP->R));
14a8d3b578SPierre Jolivet     PetscCall(VecDestroy(&lclP->lambda));
15a8d3b578SPierre Jolivet     PetscCall(VecDestroy(&lclP->lambda0));
169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->WL));
179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->W));
189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->X0));
199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->G0));
209566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GL));
219566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GAugL));
229566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->dbar));
239566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->U));
249566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->U0));
259566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->V));
269566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->V0));
279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->V1));
289566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GU));
299566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GV));
309566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GU0));
319566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GV0));
329566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GL_U));
339566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GL_V));
349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GAugL_U));
359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GAugL_V));
369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GL_U0));
379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GL_V0));
389566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GAugL_U0));
399566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->GAugL_V0));
409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->DU));
419566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->DV));
429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->WU));
439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->WV));
449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->g1));
459566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->g2));
469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->con1));
47a7e14dcfSSatish Balay 
489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->r));
499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->s));
50a7e14dcfSSatish Balay 
519566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tao->state_is));
529566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tao->design_is));
53a7e14dcfSSatish Balay 
549566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&lclP->state_scatter));
559566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&lclP->design_scatter));
56a7e14dcfSSatish Balay   }
579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lclP->R));
58a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
599566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61a7e14dcfSSatish Balay }
62a7e14dcfSSatish Balay 
TaoSetFromOptions_LCL(Tao tao,PetscOptionItems PetscOptionsObject)63ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_LCL(Tao tao, PetscOptionItems PetscOptionsObject)
64d71ae5a4SJacob Faibussowitsch {
65a7e14dcfSSatish Balay   TAO_LCL *lclP = (TAO_LCL *)tao->data;
66a7e14dcfSSatish Balay 
67a7e14dcfSSatish Balay   PetscFunctionBegin;
68d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Linearly-Constrained Augmented Lagrangian Method for PDE-constrained optimization");
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_eps1", "epsilon 1 tolerance", "", lclP->eps1, &lclP->eps1, NULL));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_eps2", "epsilon 2 tolerance", "", lclP->eps2, &lclP->eps2, NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_rho0", "init value for rho", "", lclP->rho0, &lclP->rho0, NULL));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_rhomax", "max value for rho", "", lclP->rhomax, &lclP->rhomax, NULL));
73a7e14dcfSSatish Balay   lclP->phase2_niter = 1;
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_lcl_phase2_niter", "Number of phase 2 iterations in LCL algorithm", "", lclP->phase2_niter, &lclP->phase2_niter, NULL));
75a7e14dcfSSatish Balay   lclP->verbose = PETSC_FALSE;
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_lcl_verbose", "Print verbose output", "", lclP->verbose, &lclP->verbose, NULL));
77a7e14dcfSSatish Balay   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_tola", "Tolerance for first forward solve", "", lclP->tau[0], &lclP->tau[0], NULL));
799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_tolb", "Tolerance for first adjoint solve", "", lclP->tau[1], &lclP->tau[1], NULL));
809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_tolc", "Tolerance for second forward solve", "", lclP->tau[2], &lclP->tau[2], NULL));
819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_lcl_told", "Tolerance for second adjoint solve", "", lclP->tau[3], &lclP->tau[3], NULL));
82d0609cedSBarry Smith   PetscOptionsHeadEnd();
839566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
849566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(lclP->R));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86a7e14dcfSSatish Balay }
87a7e14dcfSSatish Balay 
TaoView_LCL(Tao tao,PetscViewer viewer)88d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer)
89d71ae5a4SJacob Faibussowitsch {
903ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
91a7e14dcfSSatish Balay }
92a7e14dcfSSatish Balay 
TaoSetup_LCL(Tao tao)93d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_LCL(Tao tao)
94d71ae5a4SJacob Faibussowitsch {
95a7e14dcfSSatish Balay   TAO_LCL *lclP = (TAO_LCL *)tao->data;
96a7e14dcfSSatish Balay   PetscInt lo, hi, nlocalstate, nlocaldesign;
97a7e14dcfSSatish Balay   IS       is_state, is_design;
98f06e3bfaSBarry Smith 
99a7e14dcfSSatish Balay   PetscFunctionBegin;
1003c859ba3SBarry Smith   PetscCheck(tao->state_is, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "LCL Solver requires an initial state index set -- use TaoSetStateIS()");
1019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->gradient));
1029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
1039566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &lclP->W));
1049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &lclP->X0));
1059566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &lclP->G0));
1069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &lclP->GL));
1079566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &lclP->GAugL));
108a7e14dcfSSatish Balay 
109a8d3b578SPierre Jolivet   PetscCall(VecDuplicate(tao->constraints, &lclP->lambda));
1109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->constraints, &lclP->WL));
111a8d3b578SPierre Jolivet   PetscCall(VecDuplicate(tao->constraints, &lclP->lambda0));
1129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->constraints, &lclP->con1));
113a7e14dcfSSatish Balay 
114a8d3b578SPierre Jolivet   PetscCall(VecSet(lclP->lambda, 0.0));
115a7e14dcfSSatish Balay 
1169566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &lclP->n));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->constraints, &lclP->m));
118a7e14dcfSSatish Balay 
1199566063dSJacob Faibussowitsch   PetscCall(VecCreate(((PetscObject)tao)->comm, &lclP->U));
1209566063dSJacob Faibussowitsch   PetscCall(VecCreate(((PetscObject)tao)->comm, &lclP->V));
1219566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(tao->state_is, &nlocalstate));
1229566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(tao->design_is, &nlocaldesign));
1239566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(lclP->U, nlocalstate, lclP->m));
1249566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(lclP->V, nlocaldesign, lclP->n - lclP->m));
125f4f49eeaSPierre Jolivet   PetscCall(VecSetType(lclP->U, ((PetscObject)tao->solution)->type_name));
126f4f49eeaSPierre Jolivet   PetscCall(VecSetType(lclP->V, ((PetscObject)tao->solution)->type_name));
1279566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(lclP->U));
1289566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(lclP->V));
1299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->DU));
1309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->U0));
1319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->GU));
1329566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->GU0));
1339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->GAugL_U));
1349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->GL_U));
1359566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->GAugL_U0));
1369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->GL_U0));
1379566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->WU));
1389566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->U, &lclP->r));
1399566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->V0));
1409566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->V1));
1419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->DV));
1429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->s));
1439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->GV));
1449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->GV0));
1459566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->dbar));
1469566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->GAugL_V));
1479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->GL_V));
1489566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->GAugL_V0));
1499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->GL_V0));
1509566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->WV));
1519566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->g1));
1529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(lclP->V, &lclP->g2));
153a7e14dcfSSatish Balay 
154a7e14dcfSSatish Balay   /* create scatters for state, design subvecs */
1559566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(lclP->U, &lo, &hi));
1569566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(((PetscObject)lclP->U)->comm, hi - lo, lo, 1, &is_state));
1579566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(lclP->V, &lo, &hi));
158a7e14dcfSSatish Balay   if (0) {
159a7e14dcfSSatish Balay     PetscInt sizeU, sizeV;
1609566063dSJacob Faibussowitsch     PetscCall(VecGetSize(lclP->U, &sizeU));
1619566063dSJacob Faibussowitsch     PetscCall(VecGetSize(lclP->V, &sizeV));
16263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "size(U)=%" PetscInt_FMT ", size(V)=%" PetscInt_FMT "\n", sizeU, sizeV));
163a7e14dcfSSatish Balay   }
1649566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(((PetscObject)lclP->V)->comm, hi - lo, lo, 1, &is_design));
1659566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tao->solution, tao->state_is, lclP->U, is_state, &lclP->state_scatter));
1669566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tao->solution, tao->design_is, lclP->V, is_design, &lclP->design_scatter));
1679566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_state));
1689566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_design));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170a7e14dcfSSatish Balay }
171a7e14dcfSSatish Balay 
TaoSolve_LCL(Tao tao)172d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_LCL(Tao tao)
173d71ae5a4SJacob Faibussowitsch {
174a7e14dcfSSatish Balay   TAO_LCL                     *lclP = (TAO_LCL *)tao->data;
1758931d482SJason Sarich   PetscInt                     phase2_iter, nlocal, its;
176e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
177a7e14dcfSSatish Balay   PetscReal                    step      = 1.0, f, descent, aldescent;
178a7e14dcfSSatish Balay   PetscReal                    cnorm, mnorm;
179a7e14dcfSSatish Balay   PetscReal                    adec, r2, rGL_U, rWU;
180a7e14dcfSSatish Balay   PetscBool                    set, pset, flag, pflag, symmetric;
181a7e14dcfSSatish Balay 
182f06e3bfaSBarry Smith   PetscFunctionBegin;
183a7e14dcfSSatish Balay   lclP->rho = lclP->rho0;
1849566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(lclP->U, &nlocal));
1859566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(lclP->V, &nlocal));
1869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(lclP->R, nlocal, nlocal, lclP->n - lclP->m, lclP->n - lclP->m));
1879566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(lclP->R, lclP->V, lclP->V));
188a7e14dcfSSatish Balay   lclP->recompute_jacobian_flag = PETSC_TRUE;
189a7e14dcfSSatish Balay 
190a7e14dcfSSatish Balay   /* Scatter to U,V */
1919566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V));
192a7e14dcfSSatish Balay 
193a7e14dcfSSatish Balay   /* Evaluate Function, Gradient, Constraints, and Jacobian */
1949566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
1959566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianState(tao, tao->solution, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv));
1969566063dSJacob Faibussowitsch   PetscCall(TaoComputeJacobianDesign(tao, tao->solution, tao->jacobian_design));
1979566063dSJacob Faibussowitsch   PetscCall(TaoComputeConstraints(tao, tao->solution, tao->constraints));
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay   /* Scatter gradient to GU,GV */
2009566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV));
201a7e14dcfSSatish Balay 
202a7e14dcfSSatish Balay   /* Evaluate Lagrangian function and gradient */
203a7e14dcfSSatish Balay   /* p0 */
204a8d3b578SPierre Jolivet   PetscCall(VecSet(lclP->lambda, 0.0)); /*  Initial guess in CG */
2059566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag));
206a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
2079566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag));
208a7e14dcfSSatish Balay   } else {
209a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
210a7e14dcfSSatish Balay   }
211f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
212f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
213a7e14dcfSSatish Balay 
214a7e14dcfSSatish Balay   lclP->solve_type = LCL_ADJOINT2;
215a7e14dcfSSatish Balay   if (tao->jacobian_state_inv) {
216a7e14dcfSSatish Balay     if (symmetric) {
217a8d3b578SPierre Jolivet       PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lambda));
2189371c9d4SSatish Balay     } else {
219a8d3b578SPierre Jolivet       PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lambda));
220a7e14dcfSSatish Balay     }
221a7e14dcfSSatish Balay   } else {
2229566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre));
223a7e14dcfSSatish Balay     if (symmetric) {
224a8d3b578SPierre Jolivet       PetscCall(KSPSolve(tao->ksp, lclP->GU, lclP->lambda));
225a7e14dcfSSatish Balay     } else {
226a8d3b578SPierre Jolivet       PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lambda));
227a7e14dcfSSatish Balay     }
2289566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp, &its));
229a7e14dcfSSatish Balay     tao->ksp_its += its;
230ae93cb3cSJason Sarich     tao->ksp_tot_its += its;
231a7e14dcfSSatish Balay   }
232a8d3b578SPierre Jolivet   PetscCall(VecCopy(lclP->lambda, lclP->lambda0));
2339566063dSJacob Faibussowitsch   PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao));
234a7e14dcfSSatish Balay 
2359566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP, lclP->GL, lclP->GL_U, lclP->GL_V));
2369566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V));
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay   /* Evaluate constraint norm */
2399566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm));
2409566063dSJacob Faibussowitsch   PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm));
241a7e14dcfSSatish Balay 
242a7e14dcfSSatish Balay   /* Monitor convergence */
243763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
2449566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, mnorm, cnorm, tao->ksp_its));
2459566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, mnorm, cnorm, step));
246dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
247a7e14dcfSSatish Balay 
248763847b4SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
249e1e80dc8SAlp Dener     /* Call general purpose update function */
250dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
251ae93cb3cSJason Sarich     tao->ksp_its = 0;
252a7e14dcfSSatish Balay     /* Compute a descent direction for the linearly constrained subproblem
253a7e14dcfSSatish Balay        minimize f(u+du, v+dv)
254a7e14dcfSSatish Balay        s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay     /* Store the points around the linearization */
2579566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->U, lclP->U0));
2589566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->V, lclP->V0));
2599566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GU, lclP->GU0));
2609566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GV, lclP->GV0));
2619566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GAugL_U, lclP->GAugL_U0));
2629566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GAugL_V, lclP->GAugL_V0));
2639566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GL_U, lclP->GL_U0));
2649566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GL_V, lclP->GL_V0));
265a7e14dcfSSatish Balay 
266a7e14dcfSSatish Balay     lclP->aug0 = lclP->aug;
267a7e14dcfSSatish Balay     lclP->lgn0 = lclP->lgn;
268a7e14dcfSSatish Balay 
269a7e14dcfSSatish Balay     /* Given the design variables, we need to project the current iterate
270a7e14dcfSSatish Balay        onto the linearized constraint.  We choose to fix the design variables
271a7e14dcfSSatish Balay        and solve the linear system for the state variables.  The resulting
272a7e14dcfSSatish Balay        point is the Newton direction */
273a7e14dcfSSatish Balay 
274a7e14dcfSSatish Balay     /* Solve r = A\con */
275a7e14dcfSSatish Balay     lclP->solve_type = LCL_FORWARD1;
2769566063dSJacob Faibussowitsch     PetscCall(VecSet(lclP->r, 0.0)); /*  Initial guess in CG */
277a7e14dcfSSatish Balay 
278a7e14dcfSSatish Balay     if (tao->jacobian_state_inv) {
2799566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r));
280a7e14dcfSSatish Balay     } else {
2819566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre));
2829566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->constraints, lclP->r));
2839566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &its));
284a7e14dcfSSatish Balay       tao->ksp_its += its;
285ae93cb3cSJason Sarich       tao->ksp_tot_its += tao->ksp_its;
286a7e14dcfSSatish Balay     }
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay     /* Set design step direction dv to zero */
2899566063dSJacob Faibussowitsch     PetscCall(VecSet(lclP->s, 0.0));
290a7e14dcfSSatish Balay 
291a7e14dcfSSatish Balay     /*
292a7e14dcfSSatish Balay        Check sufficient descent for constraint merit function .5*||con||^2
293a7e14dcfSSatish Balay        con' Ak r >= eps1 ||r||^(2+eps2)
294a7e14dcfSSatish Balay     */
295a7e14dcfSSatish Balay 
296a7e14dcfSSatish Balay     /* Compute WU= Ak' * con */
297a7e14dcfSSatish Balay     if (symmetric) {
2989566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_state, tao->constraints, lclP->WU));
299a7e14dcfSSatish Balay     } else {
3009566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(tao->jacobian_state, tao->constraints, lclP->WU));
301a7e14dcfSSatish Balay     }
302a7e14dcfSSatish Balay     /* Compute r * Ak' * con */
3039566063dSJacob Faibussowitsch     PetscCall(VecDot(lclP->r, lclP->WU, &rWU));
304a7e14dcfSSatish Balay 
305a7e14dcfSSatish Balay     /* compute ||r||^(2+eps2) */
3069566063dSJacob Faibussowitsch     PetscCall(VecNorm(lclP->r, NORM_2, &r2));
307a7e14dcfSSatish Balay     r2   = PetscPowScalar(r2, 2.0 + lclP->eps2);
308a7e14dcfSSatish Balay     adec = lclP->eps1 * r2;
309a7e14dcfSSatish Balay 
310a7e14dcfSSatish Balay     if (rWU < adec) {
3119566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Newton direction not descent for constraint, feasibility phase required\n"));
31248a46eb9SPierre Jolivet       if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Newton direction not descent for constraint: %g -- using steepest descent\n", (double)descent));
313a7e14dcfSSatish Balay 
3149566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Using steepest descent direction instead.\n"));
3159566063dSJacob Faibussowitsch       PetscCall(VecSet(lclP->r, 0.0));
3169566063dSJacob Faibussowitsch       PetscCall(VecAXPY(lclP->r, -1.0, lclP->WU));
3179566063dSJacob Faibussowitsch       PetscCall(VecDot(lclP->r, lclP->r, &rWU));
3189566063dSJacob Faibussowitsch       PetscCall(VecNorm(lclP->r, NORM_2, &r2));
319a7e14dcfSSatish Balay       r2 = PetscPowScalar(r2, 2.0 + lclP->eps2);
3209566063dSJacob Faibussowitsch       PetscCall(VecDot(lclP->r, lclP->GAugL_U, &descent));
321a7e14dcfSSatish Balay       adec = lclP->eps1 * r2;
322a7e14dcfSSatish Balay     }
323a7e14dcfSSatish Balay 
324a7e14dcfSSatish Balay     /*
325a7e14dcfSSatish Balay        Check descent for aug. lagrangian
326a7e14dcfSSatish Balay        r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2)
327a7e14dcfSSatish Balay           GL_U = GUk - Ak'*yk
328a7e14dcfSSatish Balay           WU   = Ak'*con
329a7e14dcfSSatish Balay                                          adec=eps1||r||^(2+eps2)
330a7e14dcfSSatish Balay 
331a7e14dcfSSatish Balay        ==>
332a7e14dcfSSatish Balay        Check r'GL_U - rho*r'WU <= adec
333a7e14dcfSSatish Balay     */
334a7e14dcfSSatish Balay 
3359566063dSJacob Faibussowitsch     PetscCall(VecDot(lclP->r, lclP->GL_U, &rGL_U));
336a7e14dcfSSatish Balay     aldescent = rGL_U - lclP->rho * rWU;
337a7e14dcfSSatish Balay     if (aldescent > -adec) {
33848a46eb9SPierre Jolivet       if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Newton direction not descent for augmented Lagrangian: %g", (double)aldescent));
3399566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Newton direction not descent for augmented Lagrangian: %g\n", (double)aldescent));
340a7e14dcfSSatish Balay       lclP->rho = (rGL_U - adec) / rWU;
341*966bd95aSPierre Jolivet       PetscCheck(lclP->rho <= lclP->rhomax, PETSC_COMM_WORLD, PETSC_ERR_SUP, "rho=%g > rhomax, case not implemented.  Increase rhomax (-tao_lcl_rhomax)", (double)lclP->rho);
34248a46eb9SPierre Jolivet       if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Increasing penalty parameter to %g\n", (double)lclP->rho));
3439566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "  Increasing penalty parameter to %g\n", (double)lclP->rho));
344a7e14dcfSSatish Balay     }
345a7e14dcfSSatish Balay 
3469566063dSJacob Faibussowitsch     PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao));
3479566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, lclP->GL, lclP->GL_U, lclP->GL_V));
3489566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V));
349a7e14dcfSSatish Balay 
350a7e14dcfSSatish Balay     /* We now minimize the augmented Lagrangian along the Newton direction */
3519566063dSJacob Faibussowitsch     PetscCall(VecScale(lclP->r, -1.0));
3529566063dSJacob Faibussowitsch     PetscCall(LCLGather(lclP, lclP->r, lclP->s, tao->stepdirection));
3539566063dSJacob Faibussowitsch     PetscCall(VecScale(lclP->r, -1.0));
3549566063dSJacob Faibussowitsch     PetscCall(LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL));
3559566063dSJacob Faibussowitsch     PetscCall(LCLGather(lclP, lclP->U0, lclP->V0, lclP->X0));
356a7e14dcfSSatish Balay 
357a7e14dcfSSatish Balay     lclP->recompute_jacobian_flag = PETSC_TRUE;
358a7e14dcfSSatish Balay 
3599566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
3609566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
3619566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
3629566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason));
36348a46eb9SPierre Jolivet     if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steplength = %g\n", (double)step));
364a7e14dcfSSatish Balay 
3659566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V));
3669566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
3679566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV));
368a7e14dcfSSatish Balay 
3699566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V));
370a7e14dcfSSatish Balay 
371a7e14dcfSSatish Balay     /* Check convergence */
3729566063dSJacob Faibussowitsch     PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm));
3739566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm));
3749566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, mnorm, cnorm, tao->ksp_its));
3759566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, mnorm, cnorm, step));
376dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
377ad540459SPierre Jolivet     if (tao->reason != TAO_CONTINUE_ITERATING) break;
378a7e14dcfSSatish Balay 
379a7e14dcfSSatish Balay     /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */
380a7e14dcfSSatish Balay     for (phase2_iter = 0; phase2_iter < lclP->phase2_niter; phase2_iter++) {
381a7e14dcfSSatish Balay       /* We now minimize the objective function starting from the fraction of
382a7e14dcfSSatish Balay          the Newton point accepted by applying one step of a reduced-space
383a7e14dcfSSatish Balay          method.  The optimization problem is:
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay          minimize f(u+du, v+dv)
386a7e14dcfSSatish Balay          s. t.    A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0)
387a7e14dcfSSatish Balay 
388a7e14dcfSSatish Balay          In particular, we have that
389a7e14dcfSSatish Balay          du = -inv(A)*(Bdv + alpha g) */
390a7e14dcfSSatish Balay 
3919566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianState(tao, lclP->X0, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv));
3929566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianDesign(tao, lclP->X0, tao->jacobian_design));
393a7e14dcfSSatish Balay 
394a7e14dcfSSatish Balay       /* Store V and constraints */
3959566063dSJacob Faibussowitsch       PetscCall(VecCopy(lclP->V, lclP->V1));
3969566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->constraints, lclP->con1));
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay       /* Compute multipliers */
399a7e14dcfSSatish Balay       /* p1 */
400a8d3b578SPierre Jolivet       PetscCall(VecSet(lclP->lambda, 0.0)); /*  Initial guess in CG */
401a7e14dcfSSatish Balay       lclP->solve_type = LCL_ADJOINT1;
4029566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag));
403a7e14dcfSSatish Balay       if (tao->jacobian_state_pre) {
4049566063dSJacob Faibussowitsch         PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag));
405a7e14dcfSSatish Balay       } else {
406a7e14dcfSSatish Balay         pset = pflag = PETSC_TRUE;
407a7e14dcfSSatish Balay       }
408f06e3bfaSBarry Smith       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
409f06e3bfaSBarry Smith       else symmetric = PETSC_FALSE;
410a7e14dcfSSatish Balay 
411a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
412a7e14dcfSSatish Balay         if (symmetric) {
413a8d3b578SPierre Jolivet           PetscCall(MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lambda));
414a7e14dcfSSatish Balay         } else {
415a8d3b578SPierre Jolivet           PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lambda));
416a7e14dcfSSatish Balay         }
417a7e14dcfSSatish Balay       } else {
418a7e14dcfSSatish Balay         if (symmetric) {
419a8d3b578SPierre Jolivet           PetscCall(KSPSolve(tao->ksp, lclP->GAugL_U, lclP->lambda));
420a7e14dcfSSatish Balay         } else {
421a8d3b578SPierre Jolivet           PetscCall(KSPSolveTranspose(tao->ksp, lclP->GAugL_U, lclP->lambda));
422a7e14dcfSSatish Balay         }
4239566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp, &its));
424a7e14dcfSSatish Balay         tao->ksp_its += its;
425ae93cb3cSJason Sarich         tao->ksp_tot_its += its;
426a7e14dcfSSatish Balay       }
427a8d3b578SPierre Jolivet       PetscCall(MatMultTranspose(tao->jacobian_design, lclP->lambda, lclP->g1));
4289566063dSJacob Faibussowitsch       PetscCall(VecAXPY(lclP->g1, -1.0, lclP->GAugL_V));
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay       /* Compute the limited-memory quasi-newton direction */
4318931d482SJason Sarich       if (tao->niter > 0) {
4329566063dSJacob Faibussowitsch         PetscCall(MatSolve(lclP->R, lclP->g1, lclP->s));
4339566063dSJacob Faibussowitsch         PetscCall(VecDot(lclP->s, lclP->g1, &descent));
4340583ad50SJason Sarich         if (descent <= 0) {
43548a46eb9SPierre Jolivet           if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Reduced-space direction not descent: %g\n", (double)descent));
4369566063dSJacob Faibussowitsch           PetscCall(VecCopy(lclP->g1, lclP->s));
437a7e14dcfSSatish Balay         }
438a7e14dcfSSatish Balay       } else {
4399566063dSJacob Faibussowitsch         PetscCall(VecCopy(lclP->g1, lclP->s));
440a7e14dcfSSatish Balay       }
4419566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->g1, -1.0));
442a7e14dcfSSatish Balay 
443a7e14dcfSSatish Balay       /* Recover the full space direction */
4449566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_design, lclP->s, lclP->WU));
4459566063dSJacob Faibussowitsch       /* PetscCall(VecSet(lclP->r,0.0)); */ /*  Initial guess in CG */
446a7e14dcfSSatish Balay       lclP->solve_type = LCL_FORWARD2;
447a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
4489566063dSJacob Faibussowitsch         PetscCall(MatMult(tao->jacobian_state_inv, lclP->WU, lclP->r));
449a7e14dcfSSatish Balay       } else {
4509566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, lclP->WU, lclP->r));
4519566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp, &its));
452a7e14dcfSSatish Balay         tao->ksp_its += its;
453ae93cb3cSJason Sarich         tao->ksp_tot_its += its;
454a7e14dcfSSatish Balay       }
455a7e14dcfSSatish Balay 
456a7e14dcfSSatish Balay       /* We now minimize the augmented Lagrangian along the direction -r,s */
4579566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->r, -1.0));
4589566063dSJacob Faibussowitsch       PetscCall(LCLGather(lclP, lclP->r, lclP->s, tao->stepdirection));
4599566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->r, -1.0));
460a7e14dcfSSatish Balay       lclP->recompute_jacobian_flag = PETSC_TRUE;
461a7e14dcfSSatish Balay 
4629566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
4639566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
4649566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
4659566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason));
46648a46eb9SPierre Jolivet       if (lclP->verbose) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Reduced-space steplength =  %g\n", (double)step));
467a7e14dcfSSatish Balay 
4689566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V));
4699566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP, lclP->GL, lclP->GL_U, lclP->GL_V));
4709566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP, lclP->GAugL, lclP->GAugL_U, lclP->GAugL_V));
4719566063dSJacob Faibussowitsch       PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
4729566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV));
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay       /* Compute the reduced gradient at the new point */
475a7e14dcfSSatish Balay 
4769566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianState(tao, lclP->X0, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv));
4779566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianDesign(tao, lclP->X0, tao->jacobian_design));
478a7e14dcfSSatish Balay 
479a7e14dcfSSatish Balay       /* p2 */
480a8d3b578SPierre Jolivet       /* Compute multipliers, using lambda-rho*con as an initial guess in PCG */
481a7e14dcfSSatish Balay       if (phase2_iter == 0) {
482a8d3b578SPierre Jolivet         PetscCall(VecWAXPY(lclP->lambda, -lclP->rho, lclP->con1, lclP->lambda0));
483f06e3bfaSBarry Smith       } else {
484a8d3b578SPierre Jolivet         PetscCall(VecAXPY(lclP->lambda, -lclP->rho, tao->constraints));
485a7e14dcfSSatish Balay       }
486a7e14dcfSSatish Balay 
4879566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag));
488a7e14dcfSSatish Balay       if (tao->jacobian_state_pre) {
4899566063dSJacob Faibussowitsch         PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag));
490a7e14dcfSSatish Balay       } else {
491a7e14dcfSSatish Balay         pset = pflag = PETSC_TRUE;
492a7e14dcfSSatish Balay       }
493f06e3bfaSBarry Smith       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
494f06e3bfaSBarry Smith       else symmetric = PETSC_FALSE;
495a7e14dcfSSatish Balay 
496a7e14dcfSSatish Balay       lclP->solve_type = LCL_ADJOINT2;
497a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
498a7e14dcfSSatish Balay         if (symmetric) {
499a8d3b578SPierre Jolivet           PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lambda));
500a7e14dcfSSatish Balay         } else {
501a8d3b578SPierre Jolivet           PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lambda));
502a7e14dcfSSatish Balay         }
503a7e14dcfSSatish Balay       } else {
504a7e14dcfSSatish Balay         if (symmetric) {
505a8d3b578SPierre Jolivet           PetscCall(KSPSolve(tao->ksp, lclP->GU, lclP->lambda));
506a7e14dcfSSatish Balay         } else {
507a8d3b578SPierre Jolivet           PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lambda));
508a7e14dcfSSatish Balay         }
5099566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp, &its));
510a7e14dcfSSatish Balay         tao->ksp_its += its;
5112d9aa51bSJason Sarich         tao->ksp_tot_its += its;
512a7e14dcfSSatish Balay       }
513a7e14dcfSSatish Balay 
514a8d3b578SPierre Jolivet       PetscCall(MatMultTranspose(tao->jacobian_design, lclP->lambda, lclP->g2));
5159566063dSJacob Faibussowitsch       PetscCall(VecAXPY(lclP->g2, -1.0, lclP->GV));
516a7e14dcfSSatish Balay 
5179566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->g2, -1.0));
518a7e14dcfSSatish Balay 
519a7e14dcfSSatish Balay       /* Update the quasi-newton approximation */
5209566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lclP->R, lclP->V, lclP->g2));
521750b007cSBarry 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 */
522a7e14dcfSSatish Balay     }
523a7e14dcfSSatish Balay 
524a8d3b578SPierre Jolivet     PetscCall(VecCopy(lclP->lambda, lclP->lambda0));
525a7e14dcfSSatish Balay 
526a7e14dcfSSatish Balay     /* Evaluate Function, Gradient, Constraints, and Jacobian */
5279566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
5289566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, tao->solution, lclP->U, lclP->V));
5299566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP, tao->gradient, lclP->GU, lclP->GV));
530a7e14dcfSSatish Balay 
5319566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianState(tao, tao->solution, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv));
5329566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianDesign(tao, tao->solution, tao->jacobian_design));
5339566063dSJacob Faibussowitsch     PetscCall(TaoComputeConstraints(tao, tao->solution, tao->constraints));
534a7e14dcfSSatish Balay 
5359566063dSJacob Faibussowitsch     PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao));
536a7e14dcfSSatish Balay 
5379566063dSJacob Faibussowitsch     PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm));
538a7e14dcfSSatish Balay 
539a7e14dcfSSatish Balay     /* Evaluate constraint norm */
5409566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm));
541a7e14dcfSSatish Balay 
542a7e14dcfSSatish Balay     /* Monitor convergence */
5438931d482SJason Sarich     tao->niter++;
5449566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, mnorm, cnorm, tao->ksp_its));
5459566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, mnorm, cnorm, step));
546dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
547a7e14dcfSSatish Balay   }
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549a7e14dcfSSatish Balay }
550a7e14dcfSSatish Balay 
5511522df2eSJason Sarich /*MC
5521522df2eSJason Sarich  TAOLCL - linearly constrained lagrangian method for pde-constrained optimization
5531522df2eSJason Sarich 
5541522df2eSJason Sarich + -tao_lcl_eps1 - epsilon 1 tolerance
555d0609cedSBarry Smith . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);
556d0609cedSBarry Smith . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);
557d0609cedSBarry Smith . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);
5581522df2eSJason Sarich . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm
5591522df2eSJason Sarich . -tao_lcl_verbose - Print verbose output if True
5601522df2eSJason Sarich . -tao_lcl_tola - Tolerance for first forward solve
5611522df2eSJason Sarich . -tao_lcl_tolb - Tolerance for first adjoint solve
5621522df2eSJason Sarich . -tao_lcl_tolc - Tolerance for second forward solve
5631522df2eSJason Sarich - -tao_lcl_told - Tolerance for second adjoint solve
5641522df2eSJason Sarich 
5651eb8069cSJason Sarich   Level: beginner
5661522df2eSJason Sarich M*/
TaoCreate_LCL(Tao tao)567d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao)
568d71ae5a4SJacob Faibussowitsch {
569a7e14dcfSSatish Balay   TAO_LCL    *lclP;
5708caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
571a7e14dcfSSatish Balay 
572f06e3bfaSBarry Smith   PetscFunctionBegin;
573a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetup_LCL;
574a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_LCL;
575a7e14dcfSSatish Balay   tao->ops->view           = TaoView_LCL;
576a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LCL;
577a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_LCL;
5784dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lclP));
579a7e14dcfSSatish Balay   tao->data = (void *)lclP;
580a7e14dcfSSatish Balay 
5816552cf8aSJason Sarich   /* Override default settings (unless already changed) */
582606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
583606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 200);
584606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, catol, 1.0e-4);
585606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gttol, 1.0e-4);
586606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gttol, 1.0e-4);
587606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gttol, 1.0e-4);
588606f75f6SBarry Smith 
589a7e14dcfSSatish Balay   lclP->rho0       = 1.0e-4;
590a7e14dcfSSatish Balay   lclP->rhomax     = 1e5;
591a7e14dcfSSatish Balay   lclP->eps1       = 1.0e-8;
592a7e14dcfSSatish Balay   lclP->eps2       = 0.0;
593a7e14dcfSSatish Balay   lclP->solve_type = 2;
594a7e14dcfSSatish Balay   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
5959566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
5969566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
5979566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
5989566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
599a7e14dcfSSatish Balay 
6009566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, LCLComputeAugmentedLagrangianAndGradient, tao));
6019566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
6029566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
6039566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
6049566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
6050c51296cSAlp Dener 
6069566063dSJacob Faibussowitsch   PetscCall(MatCreate(((PetscObject)tao)->comm, &lclP->R));
6079566063dSJacob Faibussowitsch   PetscCall(MatSetType(lclP->R, MATLMVMBFGS));
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
609a7e14dcfSSatish Balay }
610a7e14dcfSSatish Balay 
LCLComputeLagrangianAndGradient(TaoLineSearch ls,Vec X,PetscReal * f,Vec G,void * ptr)611d71ae5a4SJacob Faibussowitsch static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
612d71ae5a4SJacob Faibussowitsch {
613441846f8SBarry Smith   Tao       tao  = (Tao)ptr;
614a7e14dcfSSatish Balay   TAO_LCL  *lclP = (TAO_LCL *)tao->data;
615a7e14dcfSSatish Balay   PetscBool set, pset, flag, pflag, symmetric;
616a7e14dcfSSatish Balay   PetscReal cdotl;
617a7e14dcfSSatish Balay 
618a7e14dcfSSatish Balay   PetscFunctionBegin;
6199566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, X, f, G));
6209566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP, G, lclP->GU, lclP->GV));
621a7e14dcfSSatish Balay   if (lclP->recompute_jacobian_flag) {
6229566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianState(tao, X, tao->jacobian_state, tao->jacobian_state_pre, tao->jacobian_state_inv));
6239566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianDesign(tao, X, tao->jacobian_design));
624a7e14dcfSSatish Balay   }
6259566063dSJacob Faibussowitsch   PetscCall(TaoComputeConstraints(tao, X, tao->constraints));
6269566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag));
627a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
6289566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag));
629a7e14dcfSSatish Balay   } else {
630a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
631a7e14dcfSSatish Balay   }
632f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
633f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
634a7e14dcfSSatish Balay 
635a8d3b578SPierre Jolivet   PetscCall(VecDot(lclP->lambda0, tao->constraints, &cdotl));
636a7e14dcfSSatish Balay   lclP->lgn = *f - cdotl;
637a7e14dcfSSatish Balay 
638a8d3b578SPierre Jolivet   /* Gradient of Lagrangian GL = G - J' * lambda */
639a7e14dcfSSatish Balay   /*      WU = A' * WL
640a7e14dcfSSatish Balay           WV = B' * WL */
641a7e14dcfSSatish Balay   if (symmetric) {
642a8d3b578SPierre Jolivet     PetscCall(MatMult(tao->jacobian_state, lclP->lambda0, lclP->GL_U));
643a7e14dcfSSatish Balay   } else {
644a8d3b578SPierre Jolivet     PetscCall(MatMultTranspose(tao->jacobian_state, lclP->lambda0, lclP->GL_U));
645a7e14dcfSSatish Balay   }
646a8d3b578SPierre Jolivet   PetscCall(MatMultTranspose(tao->jacobian_design, lclP->lambda0, lclP->GL_V));
6479566063dSJacob Faibussowitsch   PetscCall(VecScale(lclP->GL_U, -1.0));
6489566063dSJacob Faibussowitsch   PetscCall(VecScale(lclP->GL_V, -1.0));
6499566063dSJacob Faibussowitsch   PetscCall(VecAXPY(lclP->GL_U, 1.0, lclP->GU));
6509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(lclP->GL_V, 1.0, lclP->GV));
6519566063dSJacob Faibussowitsch   PetscCall(LCLGather(lclP, lclP->GL_U, lclP->GL_V, lclP->GL));
652a7e14dcfSSatish Balay 
653a7e14dcfSSatish Balay   f[0] = lclP->lgn;
6549566063dSJacob Faibussowitsch   PetscCall(VecCopy(lclP->GL, G));
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
656a7e14dcfSSatish Balay }
657a7e14dcfSSatish Balay 
LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls,Vec X,PetscReal * f,Vec G,void * ptr)658d71ae5a4SJacob Faibussowitsch static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
659d71ae5a4SJacob Faibussowitsch {
660441846f8SBarry Smith   Tao       tao  = (Tao)ptr;
661a7e14dcfSSatish Balay   TAO_LCL  *lclP = (TAO_LCL *)tao->data;
662a7e14dcfSSatish Balay   PetscReal con2;
663a7e14dcfSSatish Balay   PetscBool flag, pflag, set, pset, symmetric;
664a7e14dcfSSatish Balay 
665a7e14dcfSSatish Balay   PetscFunctionBegin;
6669566063dSJacob Faibussowitsch   PetscCall(LCLComputeLagrangianAndGradient(tao->linesearch, X, f, G, tao));
6679566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP, G, lclP->GL_U, lclP->GL_V));
6689566063dSJacob Faibussowitsch   PetscCall(VecDot(tao->constraints, tao->constraints, &con2));
669a7e14dcfSSatish Balay   lclP->aug = lclP->lgn + 0.5 * lclP->rho * con2;
670a7e14dcfSSatish Balay 
671a7e14dcfSSatish Balay   /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */
672a7e14dcfSSatish Balay   /*      WU = A' * c
673a7e14dcfSSatish Balay           WV = B' * c */
6749566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(tao->jacobian_state, &set, &flag));
675a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
6769566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre, &pset, &pflag));
677a7e14dcfSSatish Balay   } else {
678a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
679a7e14dcfSSatish Balay   }
680f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
681f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
682a7e14dcfSSatish Balay 
683a7e14dcfSSatish Balay   if (symmetric) {
6849566063dSJacob Faibussowitsch     PetscCall(MatMult(tao->jacobian_state, tao->constraints, lclP->GAugL_U));
685a7e14dcfSSatish Balay   } else {
6869566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(tao->jacobian_state, tao->constraints, lclP->GAugL_U));
687a7e14dcfSSatish Balay   }
688a7e14dcfSSatish Balay 
6899566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_design, tao->constraints, lclP->GAugL_V));
6909566063dSJacob Faibussowitsch   PetscCall(VecAYPX(lclP->GAugL_U, lclP->rho, lclP->GL_U));
6919566063dSJacob Faibussowitsch   PetscCall(VecAYPX(lclP->GAugL_V, lclP->rho, lclP->GL_V));
6929566063dSJacob Faibussowitsch   PetscCall(LCLGather(lclP, lclP->GAugL_U, lclP->GAugL_V, lclP->GAugL));
693a7e14dcfSSatish Balay 
694a7e14dcfSSatish Balay   f[0] = lclP->aug;
6959566063dSJacob Faibussowitsch   PetscCall(VecCopy(lclP->GAugL, G));
6963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
697a7e14dcfSSatish Balay }
698a7e14dcfSSatish Balay 
LCLGather(TAO_LCL * lclP,Vec u,Vec v,Vec x)69966976f2fSJacob Faibussowitsch static PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x)
700d71ae5a4SJacob Faibussowitsch {
701a7e14dcfSSatish Balay   PetscFunctionBegin;
7029566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE));
7039566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE));
7049566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE));
7059566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE));
7063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
707a7e14dcfSSatish Balay }
LCLScatter(TAO_LCL * lclP,Vec x,Vec u,Vec v)70866976f2fSJacob Faibussowitsch static PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v)
709d71ae5a4SJacob Faibussowitsch {
710a7e14dcfSSatish Balay   PetscFunctionBegin;
7119566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD));
7129566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD));
7139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD));
7149566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD));
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
716a7e14dcfSSatish Balay }
717