xref: /petsc/src/tao/pde_constrained/impls/lcl/lcl.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
7441846f8SBarry Smith static PetscErrorCode TaoDestroy_LCL(Tao tao)
8a7e14dcfSSatish Balay {
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));
149566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->lamda));
159566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lclP->lamda0));
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));
60a7e14dcfSSatish Balay   PetscFunctionReturn(0);
61a7e14dcfSSatish Balay }
62a7e14dcfSSatish Balay 
63*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_LCL(Tao tao,PetscOptionItems *PetscOptionsObject)
64a7e14dcfSSatish Balay {
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));
85a7e14dcfSSatish Balay   PetscFunctionReturn(0);
86a7e14dcfSSatish Balay }
87a7e14dcfSSatish Balay 
88441846f8SBarry Smith static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer)
89a7e14dcfSSatish Balay {
90a7e14dcfSSatish Balay   return 0;
91a7e14dcfSSatish Balay }
92a7e14dcfSSatish Balay 
93441846f8SBarry Smith static PetscErrorCode TaoSetup_LCL(Tao tao)
94a7e14dcfSSatish Balay {
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 
1099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->constraints, &lclP->lamda));
1109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->constraints, &lclP->WL));
1119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->constraints, &lclP->lamda0));
1129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->constraints, &lclP->con1));
113a7e14dcfSSatish Balay 
1149566063dSJacob Faibussowitsch   PetscCall(VecSet(lclP->lamda,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));
1259566063dSJacob Faibussowitsch   PetscCall(VecSetType(lclP->U,((PetscObject)(tao->solution))->type_name));
1269566063dSJacob Faibussowitsch   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));
169a7e14dcfSSatish Balay   PetscFunctionReturn(0);
170a7e14dcfSSatish Balay }
171a7e14dcfSSatish Balay 
172441846f8SBarry Smith static PetscErrorCode TaoSolve_LCL(Tao tao)
173a7e14dcfSSatish Balay {
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 */
2049566063dSJacob Faibussowitsch   PetscCall(VecSet(lclP->lamda,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) {
2179566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); } else {
2189566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda));
219a7e14dcfSSatish Balay     }
220a7e14dcfSSatish Balay   } else {
2219566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre));
222a7e14dcfSSatish Balay     if (symmetric) {
2239566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, lclP->GU,  lclP->lamda));
224a7e14dcfSSatish Balay     } else {
2259566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda));
226a7e14dcfSSatish Balay     }
2279566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp,&its));
228a7e14dcfSSatish Balay     tao->ksp_its+=its;
229ae93cb3cSJason Sarich     tao->ksp_tot_its+=its;
230a7e14dcfSSatish Balay   }
2319566063dSJacob Faibussowitsch   PetscCall(VecCopy(lclP->lamda,lclP->lamda0));
2329566063dSJacob Faibussowitsch   PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao));
233a7e14dcfSSatish Balay 
2349566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V));
2359566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V));
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay   /* Evaluate constraint norm */
2389566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm));
2399566063dSJacob Faibussowitsch   PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm));
240a7e14dcfSSatish Balay 
241a7e14dcfSSatish Balay   /* Monitor convergence */
242763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
2439566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its));
2449566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step));
245*dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
246a7e14dcfSSatish Balay 
247763847b4SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
248e1e80dc8SAlp Dener     /* Call general purpose update function */
249*dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao,update, tao->niter, tao->user_update);
250ae93cb3cSJason Sarich     tao->ksp_its=0;
251a7e14dcfSSatish Balay     /* Compute a descent direction for the linearly constrained subproblem
252a7e14dcfSSatish Balay        minimize f(u+du, v+dv)
253a7e14dcfSSatish Balay        s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay     /* Store the points around the linearization */
2569566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->U, lclP->U0));
2579566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->V, lclP->V0));
2589566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GU,lclP->GU0));
2599566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GV,lclP->GV0));
2609566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GAugL_U,lclP->GAugL_U0));
2619566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GAugL_V,lclP->GAugL_V0));
2629566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GL_U,lclP->GL_U0));
2639566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->GL_V,lclP->GL_V0));
264a7e14dcfSSatish Balay 
265a7e14dcfSSatish Balay     lclP->aug0 = lclP->aug;
266a7e14dcfSSatish Balay     lclP->lgn0 = lclP->lgn;
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay     /* Given the design variables, we need to project the current iterate
269a7e14dcfSSatish Balay        onto the linearized constraint.  We choose to fix the design variables
270a7e14dcfSSatish Balay        and solve the linear system for the state variables.  The resulting
271a7e14dcfSSatish Balay        point is the Newton direction */
272a7e14dcfSSatish Balay 
273a7e14dcfSSatish Balay     /* Solve r = A\con */
274a7e14dcfSSatish Balay     lclP->solve_type = LCL_FORWARD1;
2759566063dSJacob Faibussowitsch     PetscCall(VecSet(lclP->r,0.0)); /*  Initial guess in CG */
276a7e14dcfSSatish Balay 
277a7e14dcfSSatish Balay     if (tao->jacobian_state_inv) {
2789566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r));
279a7e14dcfSSatish Balay     } else {
2809566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre));
2819566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->constraints,  lclP->r));
2829566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp,&its));
283a7e14dcfSSatish Balay       tao->ksp_its+=its;
284ae93cb3cSJason Sarich       tao->ksp_tot_its+=tao->ksp_its;
285a7e14dcfSSatish Balay     }
286a7e14dcfSSatish Balay 
287a7e14dcfSSatish Balay     /* Set design step direction dv to zero */
2889566063dSJacob Faibussowitsch     PetscCall(VecSet(lclP->s, 0.0));
289a7e14dcfSSatish Balay 
290a7e14dcfSSatish Balay     /*
291a7e14dcfSSatish Balay        Check sufficient descent for constraint merit function .5*||con||^2
292a7e14dcfSSatish Balay        con' Ak r >= eps1 ||r||^(2+eps2)
293a7e14dcfSSatish Balay     */
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay     /* Compute WU= Ak' * con */
296a7e14dcfSSatish Balay     if (symmetric)  {
2979566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_state,tao->constraints,lclP->WU));
298a7e14dcfSSatish Balay     } else {
2999566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU));
300a7e14dcfSSatish Balay     }
301a7e14dcfSSatish Balay     /* Compute r * Ak' * con */
3029566063dSJacob Faibussowitsch     PetscCall(VecDot(lclP->r,lclP->WU,&rWU));
303a7e14dcfSSatish Balay 
304a7e14dcfSSatish Balay     /* compute ||r||^(2+eps2) */
3059566063dSJacob Faibussowitsch     PetscCall(VecNorm(lclP->r,NORM_2,&r2));
306a7e14dcfSSatish Balay     r2 = PetscPowScalar(r2,2.0+lclP->eps2);
307a7e14dcfSSatish Balay     adec = lclP->eps1 * r2;
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay     if (rWU < adec) {
3109566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n"));
311a7e14dcfSSatish Balay       if (lclP->verbose) {
3129566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent));
313a7e14dcfSSatish Balay       }
314a7e14dcfSSatish Balay 
3159566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"Using steepest descent direction instead.\n"));
3169566063dSJacob Faibussowitsch       PetscCall(VecSet(lclP->r,0.0));
3179566063dSJacob Faibussowitsch       PetscCall(VecAXPY(lclP->r,-1.0,lclP->WU));
3189566063dSJacob Faibussowitsch       PetscCall(VecDot(lclP->r,lclP->r,&rWU));
3199566063dSJacob Faibussowitsch       PetscCall(VecNorm(lclP->r,NORM_2,&r2));
320a7e14dcfSSatish Balay       r2 = PetscPowScalar(r2,2.0+lclP->eps2);
3219566063dSJacob Faibussowitsch       PetscCall(VecDot(lclP->r,lclP->GAugL_U,&descent));
322a7e14dcfSSatish Balay       adec = lclP->eps1 * r2;
323a7e14dcfSSatish Balay     }
324a7e14dcfSSatish Balay 
325a7e14dcfSSatish Balay     /*
326a7e14dcfSSatish Balay        Check descent for aug. lagrangian
327a7e14dcfSSatish Balay        r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2)
328a7e14dcfSSatish Balay           GL_U = GUk - Ak'*yk
329a7e14dcfSSatish Balay           WU   = Ak'*con
330a7e14dcfSSatish Balay                                          adec=eps1||r||^(2+eps2)
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay        ==>
333a7e14dcfSSatish Balay        Check r'GL_U - rho*r'WU <= adec
334a7e14dcfSSatish Balay     */
335a7e14dcfSSatish Balay 
3369566063dSJacob Faibussowitsch     PetscCall(VecDot(lclP->r,lclP->GL_U,&rGL_U));
337a7e14dcfSSatish Balay     aldescent =  rGL_U - lclP->rho*rWU;
338a7e14dcfSSatish Balay     if (aldescent > -adec) {
339a7e14dcfSSatish Balay       if (lclP->verbose) {
3409566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent));
341a7e14dcfSSatish Balay       }
3429566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"Newton direction not descent for augmented Lagrangian: %g\n",(double)aldescent));
343a7e14dcfSSatish Balay       lclP->rho =  (rGL_U - adec)/rWU;
344a7e14dcfSSatish Balay       if (lclP->rho > lclP->rhomax) {
345a7e14dcfSSatish Balay         lclP->rho = lclP->rhomax;
34698921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"rho=%g > rhomax, case not implemented.  Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho);
347a7e14dcfSSatish Balay       }
348a7e14dcfSSatish Balay       if (lclP->verbose) {
3499566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  Increasing penalty parameter to %g\n",(double)lclP->rho));
350a7e14dcfSSatish Balay       }
3519566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"  Increasing penalty parameter to %g\n",(double)lclP->rho));
352a7e14dcfSSatish Balay     }
353a7e14dcfSSatish Balay 
3549566063dSJacob Faibussowitsch     PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao));
3559566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V));
3569566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V));
357a7e14dcfSSatish Balay 
358a7e14dcfSSatish Balay     /* We now minimize the augmented Lagrangian along the Newton direction */
3599566063dSJacob Faibussowitsch     PetscCall(VecScale(lclP->r,-1.0));
3609566063dSJacob Faibussowitsch     PetscCall(LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection));
3619566063dSJacob Faibussowitsch     PetscCall(VecScale(lclP->r,-1.0));
3629566063dSJacob Faibussowitsch     PetscCall(LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL));
3639566063dSJacob Faibussowitsch     PetscCall(LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0));
364a7e14dcfSSatish Balay 
365a7e14dcfSSatish Balay     lclP->recompute_jacobian_flag = PETSC_TRUE;
366a7e14dcfSSatish Balay 
3679566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0));
3689566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
3699566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
3709566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason));
371a7e14dcfSSatish Balay     if (lclP->verbose) {
3729566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step));
373a7e14dcfSSatish Balay     }
374a7e14dcfSSatish Balay 
3759566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,tao->solution,lclP->U,lclP->V));
3769566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient));
3779566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV));
378a7e14dcfSSatish Balay 
3799566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V));
380a7e14dcfSSatish Balay 
381a7e14dcfSSatish Balay     /* Check convergence */
3829566063dSJacob Faibussowitsch     PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm));
3839566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm));
3849566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its));
3859566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step));
386*dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
387763847b4SAlp Dener     if (tao->reason != TAO_CONTINUE_ITERATING) {
388a7e14dcfSSatish Balay       break;
389a7e14dcfSSatish Balay     }
390a7e14dcfSSatish Balay 
391a7e14dcfSSatish Balay     /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */
392a7e14dcfSSatish Balay     for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++) {
393a7e14dcfSSatish Balay       /* We now minimize the objective function starting from the fraction of
394a7e14dcfSSatish Balay          the Newton point accepted by applying one step of a reduced-space
395a7e14dcfSSatish Balay          method.  The optimization problem is:
396a7e14dcfSSatish Balay 
397a7e14dcfSSatish Balay          minimize f(u+du, v+dv)
398a7e14dcfSSatish Balay          s. t.    A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0)
399a7e14dcfSSatish Balay 
400a7e14dcfSSatish Balay          In particular, we have that
401a7e14dcfSSatish Balay          du = -inv(A)*(Bdv + alpha g) */
402a7e14dcfSSatish Balay 
4039566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv));
4049566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design));
405a7e14dcfSSatish Balay 
406a7e14dcfSSatish Balay       /* Store V and constraints */
4079566063dSJacob Faibussowitsch       PetscCall(VecCopy(lclP->V, lclP->V1));
4089566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->constraints,lclP->con1));
409a7e14dcfSSatish Balay 
410a7e14dcfSSatish Balay       /* Compute multipliers */
411a7e14dcfSSatish Balay       /* p1 */
4129566063dSJacob Faibussowitsch       PetscCall(VecSet(lclP->lamda,0.0)); /*  Initial guess in CG */
413a7e14dcfSSatish Balay       lclP->solve_type = LCL_ADJOINT1;
4149566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag));
415a7e14dcfSSatish Balay       if (tao->jacobian_state_pre) {
4169566063dSJacob Faibussowitsch         PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag));
417a7e14dcfSSatish Balay       } else {
418a7e14dcfSSatish Balay         pset = pflag = PETSC_TRUE;
419a7e14dcfSSatish Balay       }
420f06e3bfaSBarry Smith       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
421f06e3bfaSBarry Smith       else symmetric = PETSC_FALSE;
422a7e14dcfSSatish Balay 
423a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
424a7e14dcfSSatish Balay         if (symmetric) {
4259566063dSJacob Faibussowitsch           PetscCall(MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda));
426a7e14dcfSSatish Balay         } else {
4279566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda));
428a7e14dcfSSatish Balay         }
429a7e14dcfSSatish Balay       } else {
430a7e14dcfSSatish Balay         if (symmetric) {
4319566063dSJacob Faibussowitsch           PetscCall(KSPSolve(tao->ksp, lclP->GAugL_U,  lclP->lamda));
432a7e14dcfSSatish Balay         } else {
4339566063dSJacob Faibussowitsch           PetscCall(KSPSolveTranspose(tao->ksp, lclP->GAugL_U,  lclP->lamda));
434a7e14dcfSSatish Balay         }
4359566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp,&its));
436a7e14dcfSSatish Balay         tao->ksp_its+=its;
437ae93cb3cSJason Sarich         tao->ksp_tot_its+=its;
438a7e14dcfSSatish Balay       }
4399566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1));
4409566063dSJacob Faibussowitsch       PetscCall(VecAXPY(lclP->g1,-1.0,lclP->GAugL_V));
441a7e14dcfSSatish Balay 
442a7e14dcfSSatish Balay       /* Compute the limited-memory quasi-newton direction */
4438931d482SJason Sarich       if (tao->niter > 0) {
4449566063dSJacob Faibussowitsch         PetscCall(MatSolve(lclP->R,lclP->g1,lclP->s));
4459566063dSJacob Faibussowitsch         PetscCall(VecDot(lclP->s,lclP->g1,&descent));
4460583ad50SJason Sarich         if (descent <= 0) {
447a7e14dcfSSatish Balay           if (lclP->verbose) {
4489566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent));
449a7e14dcfSSatish Balay           }
4509566063dSJacob Faibussowitsch           PetscCall(VecCopy(lclP->g1,lclP->s));
451a7e14dcfSSatish Balay         }
452a7e14dcfSSatish Balay       } else {
4539566063dSJacob Faibussowitsch         PetscCall(VecCopy(lclP->g1,lclP->s));
454a7e14dcfSSatish Balay       }
4559566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->g1,-1.0));
456a7e14dcfSSatish Balay 
457a7e14dcfSSatish Balay       /* Recover the full space direction */
4589566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->jacobian_design,lclP->s,lclP->WU));
4599566063dSJacob Faibussowitsch       /* PetscCall(VecSet(lclP->r,0.0)); */ /*  Initial guess in CG */
460a7e14dcfSSatish Balay       lclP->solve_type = LCL_FORWARD2;
461a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
4629566063dSJacob Faibussowitsch         PetscCall(MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r));
463a7e14dcfSSatish Balay       } else {
4649566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, lclP->WU, lclP->r));
4659566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp,&its));
466a7e14dcfSSatish Balay         tao->ksp_its += its;
467ae93cb3cSJason Sarich         tao->ksp_tot_its+=its;
468a7e14dcfSSatish Balay       }
469a7e14dcfSSatish Balay 
470a7e14dcfSSatish Balay       /* We now minimize the augmented Lagrangian along the direction -r,s */
4719566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->r, -1.0));
4729566063dSJacob Faibussowitsch       PetscCall(LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection));
4739566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->r, -1.0));
474a7e14dcfSSatish Balay       lclP->recompute_jacobian_flag = PETSC_TRUE;
475a7e14dcfSSatish Balay 
4769566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0));
4779566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT));
4789566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
4799566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason));
480a7e14dcfSSatish Balay       if (lclP->verbose) {
4819566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength =  %g\n",(double)step));
482a7e14dcfSSatish Balay       }
483a7e14dcfSSatish Balay 
4849566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP,tao->solution,lclP->U,lclP->V));
4859566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V));
4869566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V));
4879566063dSJacob Faibussowitsch       PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient));
4889566063dSJacob Faibussowitsch       PetscCall(LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV));
489a7e14dcfSSatish Balay 
490a7e14dcfSSatish Balay       /* Compute the reduced gradient at the new point */
491a7e14dcfSSatish Balay 
4929566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv));
4939566063dSJacob Faibussowitsch       PetscCall(TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design));
494a7e14dcfSSatish Balay 
495a7e14dcfSSatish Balay       /* p2 */
496a7e14dcfSSatish Balay       /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */
497a7e14dcfSSatish Balay       if (phase2_iter==0) {
4989566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0));
499f06e3bfaSBarry Smith       } else {
5009566063dSJacob Faibussowitsch         PetscCall(VecAXPY(lclP->lamda,-lclP->rho,tao->constraints));
501a7e14dcfSSatish Balay       }
502a7e14dcfSSatish Balay 
5039566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag));
504a7e14dcfSSatish Balay       if (tao->jacobian_state_pre) {
5059566063dSJacob Faibussowitsch         PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag));
506a7e14dcfSSatish Balay       } else {
507a7e14dcfSSatish Balay         pset = pflag = PETSC_TRUE;
508a7e14dcfSSatish Balay       }
509f06e3bfaSBarry Smith       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
510f06e3bfaSBarry Smith       else symmetric = PETSC_FALSE;
511a7e14dcfSSatish Balay 
512a7e14dcfSSatish Balay       lclP->solve_type = LCL_ADJOINT2;
513a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
514a7e14dcfSSatish Balay         if (symmetric) {
5159566063dSJacob Faibussowitsch           PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda));
516a7e14dcfSSatish Balay         } else {
5179566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda));
518a7e14dcfSSatish Balay         }
519a7e14dcfSSatish Balay       } else {
520a7e14dcfSSatish Balay         if (symmetric) {
5219566063dSJacob Faibussowitsch           PetscCall(KSPSolve(tao->ksp, lclP->GU,  lclP->lamda));
522a7e14dcfSSatish Balay         } else {
5239566063dSJacob Faibussowitsch           PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda));
524a7e14dcfSSatish Balay         }
5259566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp,&its));
526a7e14dcfSSatish Balay         tao->ksp_its += its;
5272d9aa51bSJason Sarich         tao->ksp_tot_its += its;
528a7e14dcfSSatish Balay       }
529a7e14dcfSSatish Balay 
5309566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2));
5319566063dSJacob Faibussowitsch       PetscCall(VecAXPY(lclP->g2,-1.0,lclP->GV));
532a7e14dcfSSatish Balay 
5339566063dSJacob Faibussowitsch       PetscCall(VecScale(lclP->g2,-1.0));
534a7e14dcfSSatish Balay 
535a7e14dcfSSatish Balay       /* Update the quasi-newton approximation */
5369566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lclP->R,lclP->V,lclP->g2));
537a7e14dcfSSatish Balay       /* 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 */
538a7e14dcfSSatish Balay 
539a7e14dcfSSatish Balay     }
540a7e14dcfSSatish Balay 
5419566063dSJacob Faibussowitsch     PetscCall(VecCopy(lclP->lamda,lclP->lamda0));
542a7e14dcfSSatish Balay 
543a7e14dcfSSatish Balay     /* Evaluate Function, Gradient, Constraints, and Jacobian */
5449566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient));
5459566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,tao->solution,lclP->U,lclP->V));
5469566063dSJacob Faibussowitsch     PetscCall(LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV));
547a7e14dcfSSatish Balay 
5489566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv));
5499566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design));
5509566063dSJacob Faibussowitsch     PetscCall(TaoComputeConstraints(tao,tao->solution, tao->constraints));
551a7e14dcfSSatish Balay 
5529566063dSJacob Faibussowitsch     PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao));
553a7e14dcfSSatish Balay 
5549566063dSJacob Faibussowitsch     PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm));
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay     /* Evaluate constraint norm */
5579566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm));
558a7e14dcfSSatish Balay 
559a7e14dcfSSatish Balay     /* Monitor convergence */
5608931d482SJason Sarich     tao->niter++;
5619566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its));
5629566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step));
563*dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
564a7e14dcfSSatish Balay   }
565a7e14dcfSSatish Balay   PetscFunctionReturn(0);
566a7e14dcfSSatish Balay }
567a7e14dcfSSatish Balay 
5681522df2eSJason Sarich /*MC
5691522df2eSJason Sarich  TAOLCL - linearly constrained lagrangian method for pde-constrained optimization
5701522df2eSJason Sarich 
5711522df2eSJason Sarich + -tao_lcl_eps1 - epsilon 1 tolerance
572d0609cedSBarry Smith . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);
573d0609cedSBarry Smith . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);
574d0609cedSBarry Smith . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);
5751522df2eSJason Sarich . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm
5761522df2eSJason Sarich . -tao_lcl_verbose - Print verbose output if True
5771522df2eSJason Sarich . -tao_lcl_tola - Tolerance for first forward solve
5781522df2eSJason Sarich . -tao_lcl_tolb - Tolerance for first adjoint solve
5791522df2eSJason Sarich . -tao_lcl_tolc - Tolerance for second forward solve
5801522df2eSJason Sarich - -tao_lcl_told - Tolerance for second adjoint solve
5811522df2eSJason Sarich 
5821eb8069cSJason Sarich   Level: beginner
5831522df2eSJason Sarich M*/
584728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao)
585a7e14dcfSSatish Balay {
586a7e14dcfSSatish Balay   TAO_LCL        *lclP;
5878caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
588a7e14dcfSSatish Balay 
589f06e3bfaSBarry Smith   PetscFunctionBegin;
590a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_LCL;
591a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_LCL;
592a7e14dcfSSatish Balay   tao->ops->view = TaoView_LCL;
593a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LCL;
594a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_LCL;
5959566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&lclP));
596a7e14dcfSSatish Balay   tao->data = (void*)lclP;
597a7e14dcfSSatish Balay 
5986552cf8aSJason Sarich   /* Override default settings (unless already changed) */
5996552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 200;
6006552cf8aSJason Sarich   if (!tao->catol_changed) tao->catol = 1.0e-4;
6016552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gttol = 1.0e-4;
6026552cf8aSJason Sarich   if (!tao->grtol_changed) tao->gttol = 1.0e-4;
6036552cf8aSJason Sarich   if (!tao->gttol_changed) tao->gttol = 1.0e-4;
604a7e14dcfSSatish Balay   lclP->rho0 = 1.0e-4;
605a7e14dcfSSatish Balay   lclP->rhomax=1e5;
606a7e14dcfSSatish Balay   lclP->eps1 = 1.0e-8;
607a7e14dcfSSatish Balay   lclP->eps2 = 0.0;
608a7e14dcfSSatish Balay   lclP->solve_type=2;
609a7e14dcfSSatish Balay   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
6109566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
6119566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
6129566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
6139566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
614a7e14dcfSSatish Balay 
6159566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao));
6169566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp));
6179566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
6189566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
6199566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
6200c51296cSAlp Dener 
6219566063dSJacob Faibussowitsch   PetscCall(MatCreate(((PetscObject)tao)->comm, &lclP->R));
6229566063dSJacob Faibussowitsch   PetscCall(MatSetType(lclP->R, MATLMVMBFGS));
623a7e14dcfSSatish Balay   PetscFunctionReturn(0);
624a7e14dcfSSatish Balay }
625a7e14dcfSSatish Balay 
626a7e14dcfSSatish Balay static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
627a7e14dcfSSatish Balay {
628441846f8SBarry Smith   Tao            tao = (Tao)ptr;
629a7e14dcfSSatish Balay   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
630a7e14dcfSSatish Balay   PetscBool      set,pset,flag,pflag,symmetric;
631a7e14dcfSSatish Balay   PetscReal      cdotl;
632a7e14dcfSSatish Balay 
633a7e14dcfSSatish Balay   PetscFunctionBegin;
6349566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao,X,f,G));
6359566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP,G,lclP->GU,lclP->GV));
636a7e14dcfSSatish Balay   if (lclP->recompute_jacobian_flag) {
6379566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv));
6389566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianDesign(tao,X,tao->jacobian_design));
639a7e14dcfSSatish Balay   }
6409566063dSJacob Faibussowitsch   PetscCall(TaoComputeConstraints(tao,X, tao->constraints));
6419566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag));
642a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
6439566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag));
644a7e14dcfSSatish Balay   } else {
645a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
646a7e14dcfSSatish Balay   }
647f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
648f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
649a7e14dcfSSatish Balay 
6509566063dSJacob Faibussowitsch   PetscCall(VecDot(lclP->lamda0, tao->constraints, &cdotl));
651a7e14dcfSSatish Balay   lclP->lgn = *f - cdotl;
652a7e14dcfSSatish Balay 
653a7e14dcfSSatish Balay   /* Gradient of Lagrangian GL = G - J' * lamda */
654a7e14dcfSSatish Balay   /*      WU = A' * WL
655a7e14dcfSSatish Balay           WV = B' * WL */
656a7e14dcfSSatish Balay   if (symmetric) {
6579566063dSJacob Faibussowitsch     PetscCall(MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U));
658a7e14dcfSSatish Balay   } else {
6599566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U));
660a7e14dcfSSatish Balay   }
6619566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V));
6629566063dSJacob Faibussowitsch   PetscCall(VecScale(lclP->GL_U,-1.0));
6639566063dSJacob Faibussowitsch   PetscCall(VecScale(lclP->GL_V,-1.0));
6649566063dSJacob Faibussowitsch   PetscCall(VecAXPY(lclP->GL_U,1.0,lclP->GU));
6659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(lclP->GL_V,1.0,lclP->GV));
6669566063dSJacob Faibussowitsch   PetscCall(LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL));
667a7e14dcfSSatish Balay 
668a7e14dcfSSatish Balay   f[0] = lclP->lgn;
6699566063dSJacob Faibussowitsch   PetscCall(VecCopy(lclP->GL,G));
670a7e14dcfSSatish Balay   PetscFunctionReturn(0);
671a7e14dcfSSatish Balay }
672a7e14dcfSSatish Balay 
673a7e14dcfSSatish Balay static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
674a7e14dcfSSatish Balay {
675441846f8SBarry Smith   Tao            tao = (Tao)ptr;
676a7e14dcfSSatish Balay   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
677a7e14dcfSSatish Balay   PetscReal      con2;
678a7e14dcfSSatish Balay   PetscBool      flag,pflag,set,pset,symmetric;
679a7e14dcfSSatish Balay 
680a7e14dcfSSatish Balay   PetscFunctionBegin;
6819566063dSJacob Faibussowitsch   PetscCall(LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao));
6829566063dSJacob Faibussowitsch   PetscCall(LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V));
6839566063dSJacob Faibussowitsch   PetscCall(VecDot(tao->constraints,tao->constraints,&con2));
684a7e14dcfSSatish Balay   lclP->aug = lclP->lgn + 0.5*lclP->rho*con2;
685a7e14dcfSSatish Balay 
686a7e14dcfSSatish Balay   /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */
687a7e14dcfSSatish Balay   /*      WU = A' * c
688a7e14dcfSSatish Balay           WV = B' * c */
6899566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag));
690a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
6919566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag));
692a7e14dcfSSatish Balay   } else {
693a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
694a7e14dcfSSatish Balay   }
695f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
696f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
697a7e14dcfSSatish Balay 
698a7e14dcfSSatish Balay   if (symmetric) {
6999566063dSJacob Faibussowitsch     PetscCall(MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U));
700a7e14dcfSSatish Balay   } else {
7019566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U));
702a7e14dcfSSatish Balay   }
703a7e14dcfSSatish Balay 
7049566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V));
7059566063dSJacob Faibussowitsch   PetscCall(VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U));
7069566063dSJacob Faibussowitsch   PetscCall(VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V));
7079566063dSJacob Faibussowitsch   PetscCall(LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL));
708a7e14dcfSSatish Balay 
709a7e14dcfSSatish Balay   f[0] = lclP->aug;
7109566063dSJacob Faibussowitsch   PetscCall(VecCopy(lclP->GAugL,G));
711a7e14dcfSSatish Balay   PetscFunctionReturn(0);
712a7e14dcfSSatish Balay }
713a7e14dcfSSatish Balay 
714a7e14dcfSSatish Balay PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x)
715a7e14dcfSSatish Balay {
716a7e14dcfSSatish Balay   PetscFunctionBegin;
7179566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE));
7189566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE));
7199566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE));
7209566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE));
721a7e14dcfSSatish Balay   PetscFunctionReturn(0);
722a7e14dcfSSatish Balay 
723a7e14dcfSSatish Balay }
724a7e14dcfSSatish Balay PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v)
725a7e14dcfSSatish Balay {
726a7e14dcfSSatish Balay   PetscFunctionBegin;
7279566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD));
7289566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD));
7299566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD));
7309566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD));
731a7e14dcfSSatish Balay   PetscFunctionReturn(0);
732a7e14dcfSSatish Balay 
733a7e14dcfSSatish Balay }
734