xref: /petsc/src/tao/pde_constrained/impls/lcl/lcl.c (revision 3c859ba3a04a72e8efcb87bc7ffc046a6cbab413)
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;
10a7e14dcfSSatish Balay   PetscErrorCode ierr;
11f06e3bfaSBarry Smith 
12a7e14dcfSSatish Balay   PetscFunctionBegin;
13a7e14dcfSSatish Balay   if (tao->setupcalled) {
14a7e14dcfSSatish Balay     ierr = MatDestroy(&lclP->R);CHKERRQ(ierr);
15a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->lamda);CHKERRQ(ierr);
16a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->lamda0);CHKERRQ(ierr);
17a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->WL);CHKERRQ(ierr);
18a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->W);CHKERRQ(ierr);
19a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->X0);CHKERRQ(ierr);
20a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->G0);CHKERRQ(ierr);
21a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GL);CHKERRQ(ierr);
22a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GAugL);CHKERRQ(ierr);
23a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->dbar);CHKERRQ(ierr);
24a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->U);CHKERRQ(ierr);
25a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->U0);CHKERRQ(ierr);
26a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->V);CHKERRQ(ierr);
27a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->V0);CHKERRQ(ierr);
28a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->V1);CHKERRQ(ierr);
29a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GU);CHKERRQ(ierr);
30a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GV);CHKERRQ(ierr);
31a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GU0);CHKERRQ(ierr);
32a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GV0);CHKERRQ(ierr);
33a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GL_U);CHKERRQ(ierr);
34a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GL_V);CHKERRQ(ierr);
35a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GAugL_U);CHKERRQ(ierr);
36a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GAugL_V);CHKERRQ(ierr);
37a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GL_U0);CHKERRQ(ierr);
38a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GL_V0);CHKERRQ(ierr);
39a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GAugL_U0);CHKERRQ(ierr);
40a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->GAugL_V0);CHKERRQ(ierr);
41a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->DU);CHKERRQ(ierr);
42a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->DV);CHKERRQ(ierr);
43a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->WU);CHKERRQ(ierr);
44a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->WV);CHKERRQ(ierr);
45a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->g1);CHKERRQ(ierr);
46a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->g2);CHKERRQ(ierr);
47a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->con1);CHKERRQ(ierr);
48a7e14dcfSSatish Balay 
49a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->r);CHKERRQ(ierr);
50a7e14dcfSSatish Balay     ierr = VecDestroy(&lclP->s);CHKERRQ(ierr);
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay     ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr);
53a7e14dcfSSatish Balay     ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr);
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay     ierr = VecScatterDestroy(&lclP->state_scatter);CHKERRQ(ierr);
56a7e14dcfSSatish Balay     ierr = VecScatterDestroy(&lclP->design_scatter);CHKERRQ(ierr);
57a7e14dcfSSatish Balay   }
580c51296cSAlp Dener   ierr = MatDestroy(&lclP->R);CHKERRQ(ierr);
59302440fdSBarry Smith   ierr = PetscFree(tao->data);CHKERRQ(ierr);
60a7e14dcfSSatish Balay   PetscFunctionReturn(0);
61a7e14dcfSSatish Balay }
62a7e14dcfSSatish Balay 
634416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LCL(PetscOptionItems *PetscOptionsObject,Tao tao)
64a7e14dcfSSatish Balay {
65a7e14dcfSSatish Balay   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
66a7e14dcfSSatish Balay   PetscErrorCode ierr;
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay   PetscFunctionBegin;
691a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Linearly-Constrained Augmented Lagrangian Method for PDE-constrained optimization");CHKERRQ(ierr);
7094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_eps1","epsilon 1 tolerance","",lclP->eps1,&lclP->eps1,NULL);CHKERRQ(ierr);
7194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr);
7294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr);
7394ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr);
74a7e14dcfSSatish Balay   lclP->phase2_niter = 1;
7594ae4db5SBarry Smith   ierr = PetscOptionsInt("-tao_lcl_phase2_niter","Number of phase 2 iterations in LCL algorithm","",lclP->phase2_niter,&lclP->phase2_niter,NULL);CHKERRQ(ierr);
76a7e14dcfSSatish Balay   lclP->verbose = PETSC_FALSE;
7794ae4db5SBarry Smith   ierr = PetscOptionsBool("-tao_lcl_verbose","Print verbose output","",lclP->verbose,&lclP->verbose,NULL);CHKERRQ(ierr);
78a7e14dcfSSatish Balay   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
7994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_tola","Tolerance for first forward solve","",lclP->tau[0],&lclP->tau[0],NULL);CHKERRQ(ierr);
8094ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_tolb","Tolerance for first adjoint solve","",lclP->tau[1],&lclP->tau[1],NULL);CHKERRQ(ierr);
8194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_tolc","Tolerance for second forward solve","",lclP->tau[2],&lclP->tau[2],NULL);CHKERRQ(ierr);
8294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_lcl_told","Tolerance for second adjoint solve","",lclP->tau[3],&lclP->tau[3],NULL);CHKERRQ(ierr);
83a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
84a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
850c51296cSAlp Dener   ierr = MatSetFromOptions(lclP->R);CHKERRQ(ierr);
86a7e14dcfSSatish Balay   PetscFunctionReturn(0);
87a7e14dcfSSatish Balay }
88a7e14dcfSSatish Balay 
89441846f8SBarry Smith static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer)
90a7e14dcfSSatish Balay {
91a7e14dcfSSatish Balay   return 0;
92a7e14dcfSSatish Balay }
93a7e14dcfSSatish Balay 
94441846f8SBarry Smith static PetscErrorCode TaoSetup_LCL(Tao tao)
95a7e14dcfSSatish Balay {
96a7e14dcfSSatish Balay   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
97a7e14dcfSSatish Balay   PetscInt       lo, hi, nlocalstate, nlocaldesign;
98a7e14dcfSSatish Balay   PetscErrorCode ierr;
99a7e14dcfSSatish Balay   IS             is_state, is_design;
100f06e3bfaSBarry Smith 
101a7e14dcfSSatish Balay   PetscFunctionBegin;
102*3c859ba3SBarry Smith   PetscCheck(tao->state_is,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"LCL Solver requires an initial state index set -- use TaoSetStateIS()");
103a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);
104a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
105a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &lclP->W);CHKERRQ(ierr);
106a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &lclP->X0);CHKERRQ(ierr);
107a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &lclP->G0);CHKERRQ(ierr);
108a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &lclP->GL);CHKERRQ(ierr);
109a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &lclP->GAugL);CHKERRQ(ierr);
110a7e14dcfSSatish Balay 
111a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->constraints, &lclP->lamda);CHKERRQ(ierr);
112a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->constraints, &lclP->WL);CHKERRQ(ierr);
113a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->constraints, &lclP->lamda0);CHKERRQ(ierr);
114a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->constraints, &lclP->con1);CHKERRQ(ierr);
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay   ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr);
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution, &lclP->n);CHKERRQ(ierr);
119a7e14dcfSSatish Balay   ierr = VecGetSize(tao->constraints, &lclP->m);CHKERRQ(ierr);
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   ierr = VecCreate(((PetscObject)tao)->comm,&lclP->U);CHKERRQ(ierr);
122a7e14dcfSSatish Balay   ierr = VecCreate(((PetscObject)tao)->comm,&lclP->V);CHKERRQ(ierr);
123a7e14dcfSSatish Balay   ierr = ISGetLocalSize(tao->state_is,&nlocalstate);CHKERRQ(ierr);
124a7e14dcfSSatish Balay   ierr = ISGetLocalSize(tao->design_is,&nlocaldesign);CHKERRQ(ierr);
125a7e14dcfSSatish Balay   ierr = VecSetSizes(lclP->U,nlocalstate,lclP->m);CHKERRQ(ierr);
126a7e14dcfSSatish Balay   ierr = VecSetSizes(lclP->V,nlocaldesign,lclP->n-lclP->m);CHKERRQ(ierr);
127a7e14dcfSSatish Balay   ierr = VecSetType(lclP->U,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr);
128a7e14dcfSSatish Balay   ierr = VecSetType(lclP->V,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr);
129a7e14dcfSSatish Balay   ierr = VecSetFromOptions(lclP->U);CHKERRQ(ierr);
130a7e14dcfSSatish Balay   ierr = VecSetFromOptions(lclP->V);CHKERRQ(ierr);
131a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->DU);CHKERRQ(ierr);
132a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->U0);CHKERRQ(ierr);
133a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->GU);CHKERRQ(ierr);
134a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->GU0);CHKERRQ(ierr);
135a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->GAugL_U);CHKERRQ(ierr);
136a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->GL_U);CHKERRQ(ierr);
137a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->GAugL_U0);CHKERRQ(ierr);
138a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->GL_U0);CHKERRQ(ierr);
139a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->WU);CHKERRQ(ierr);
140a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->U,&lclP->r);CHKERRQ(ierr);
141a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->V0);CHKERRQ(ierr);
142a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->V1);CHKERRQ(ierr);
143a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->DV);CHKERRQ(ierr);
144a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->s);CHKERRQ(ierr);
145a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->GV);CHKERRQ(ierr);
146a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->GV0);CHKERRQ(ierr);
147a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->dbar);CHKERRQ(ierr);
148a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->GAugL_V);CHKERRQ(ierr);
149a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->GL_V);CHKERRQ(ierr);
150a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->GAugL_V0);CHKERRQ(ierr);
151a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->GL_V0);CHKERRQ(ierr);
152a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->WV);CHKERRQ(ierr);
153a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->g1);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   ierr = VecDuplicate(lclP->V,&lclP->g2);CHKERRQ(ierr);
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay   /* create scatters for state, design subvecs */
157a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(lclP->U,&lo,&hi);CHKERRQ(ierr);
158a7e14dcfSSatish Balay   ierr = ISCreateStride(((PetscObject)lclP->U)->comm,hi-lo,lo,1,&is_state);CHKERRQ(ierr);
159a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(lclP->V,&lo,&hi);CHKERRQ(ierr);
160a7e14dcfSSatish Balay   if (0) {
161a7e14dcfSSatish Balay     PetscInt sizeU,sizeV;
162928bb9adSStefano Zampini     ierr = VecGetSize(lclP->U,&sizeU);CHKERRQ(ierr);
163928bb9adSStefano Zampini     ierr = VecGetSize(lclP->V,&sizeV);CHKERRQ(ierr);
164928bb9adSStefano Zampini     ierr = PetscPrintf(PETSC_COMM_WORLD,"size(U)=%D, size(V)=%D\n",sizeU,sizeV);CHKERRQ(ierr);
165a7e14dcfSSatish Balay   }
166a7e14dcfSSatish Balay   ierr = ISCreateStride(((PetscObject)lclP->V)->comm,hi-lo,lo,1,&is_design);CHKERRQ(ierr);
1679448b7f1SJunchao Zhang   ierr = VecScatterCreate(tao->solution,tao->state_is,lclP->U,is_state,&lclP->state_scatter);CHKERRQ(ierr);
1689448b7f1SJunchao Zhang   ierr = VecScatterCreate(tao->solution,tao->design_is,lclP->V,is_design,&lclP->design_scatter);CHKERRQ(ierr);
169a7e14dcfSSatish Balay   ierr = ISDestroy(&is_state);CHKERRQ(ierr);
170a7e14dcfSSatish Balay   ierr = ISDestroy(&is_design);CHKERRQ(ierr);
171a7e14dcfSSatish Balay   PetscFunctionReturn(0);
172a7e14dcfSSatish Balay }
173a7e14dcfSSatish Balay 
174441846f8SBarry Smith static PetscErrorCode TaoSolve_LCL(Tao tao)
175a7e14dcfSSatish Balay {
176a7e14dcfSSatish Balay   TAO_LCL                      *lclP = (TAO_LCL*)tao->data;
1778931d482SJason Sarich   PetscInt                     phase2_iter,nlocal,its;
178e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
179a7e14dcfSSatish Balay   PetscReal                    step=1.0,f, descent, aldescent;
180a7e14dcfSSatish Balay   PetscReal                    cnorm, mnorm;
181a7e14dcfSSatish Balay   PetscReal                    adec,r2,rGL_U,rWU;
182a7e14dcfSSatish Balay   PetscBool                    set,pset,flag,pflag,symmetric;
183a7e14dcfSSatish Balay   PetscErrorCode               ierr;
184a7e14dcfSSatish Balay 
185f06e3bfaSBarry Smith   PetscFunctionBegin;
186a7e14dcfSSatish Balay   lclP->rho = lclP->rho0;
187a7e14dcfSSatish Balay   ierr = VecGetLocalSize(lclP->U,&nlocal);CHKERRQ(ierr);
188a7e14dcfSSatish Balay   ierr = VecGetLocalSize(lclP->V,&nlocal);CHKERRQ(ierr);
189f5a7d1c1SBarry Smith   ierr = MatSetSizes(lclP->R, nlocal, nlocal, lclP->n-lclP->m, lclP->n-lclP->m);CHKERRQ(ierr);
190cd929ea3SAlp Dener   ierr = MatLMVMAllocate(lclP->R,lclP->V,lclP->V);CHKERRQ(ierr);
191a7e14dcfSSatish Balay   lclP->recompute_jacobian_flag = PETSC_TRUE;
192a7e14dcfSSatish Balay 
193a7e14dcfSSatish Balay   /* Scatter to U,V */
194a7e14dcfSSatish Balay   ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay   /* Evaluate Function, Gradient, Constraints, and Jacobian */
197a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
198ffad9901SBarry Smith   ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
19994ab13aaSBarry Smith   ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr);
200a7e14dcfSSatish Balay   ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr);
201a7e14dcfSSatish Balay 
202a7e14dcfSSatish Balay   /* Scatter gradient to GU,GV */
203a7e14dcfSSatish Balay   ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
204a7e14dcfSSatish Balay 
205a7e14dcfSSatish Balay   /* Evaluate Lagrangian function and gradient */
206a7e14dcfSSatish Balay   /* p0 */
207a7e14dcfSSatish Balay   ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
208a7e14dcfSSatish Balay   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
209a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
210a7e14dcfSSatish Balay     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
211a7e14dcfSSatish Balay   } else {
212a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
213a7e14dcfSSatish Balay   }
214f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
215f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
216a7e14dcfSSatish Balay 
217a7e14dcfSSatish Balay   lclP->solve_type = LCL_ADJOINT2;
218a7e14dcfSSatish Balay   if (tao->jacobian_state_inv) {
219a7e14dcfSSatish Balay     if (symmetric) {
220a7e14dcfSSatish Balay       ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); } else {
221a7e14dcfSSatish Balay       ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
222a7e14dcfSSatish Balay     }
223a7e14dcfSSatish Balay   } else {
22423ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr);
225a7e14dcfSSatish Balay     if (symmetric) {
226a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
227a7e14dcfSSatish Balay     } else {
228a7e14dcfSSatish Balay       ierr = KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
229a7e14dcfSSatish Balay     }
230a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
231a7e14dcfSSatish Balay     tao->ksp_its+=its;
232ae93cb3cSJason Sarich     tao->ksp_tot_its+=its;
233a7e14dcfSSatish Balay   }
234a7e14dcfSSatish Balay   ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr);
235a7e14dcfSSatish Balay   ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay   ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
238a7e14dcfSSatish Balay   ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay   /* Evaluate constraint norm */
241a7e14dcfSSatish Balay   ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
242a7e14dcfSSatish Balay   ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay   /* Monitor convergence */
245763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
246763847b4SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
247763847b4SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
248763847b4SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
249a7e14dcfSSatish Balay 
250763847b4SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
251e1e80dc8SAlp Dener     /* Call general purpose update function */
252e1e80dc8SAlp Dener     if (tao->ops->update) {
2538fcddce6SStefano Zampini       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
254e1e80dc8SAlp Dener     }
255ae93cb3cSJason Sarich     tao->ksp_its=0;
256a7e14dcfSSatish Balay     /* Compute a descent direction for the linearly constrained subproblem
257a7e14dcfSSatish Balay        minimize f(u+du, v+dv)
258a7e14dcfSSatish Balay        s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay     /* Store the points around the linearization */
261a7e14dcfSSatish Balay     ierr = VecCopy(lclP->U, lclP->U0);CHKERRQ(ierr);
262a7e14dcfSSatish Balay     ierr = VecCopy(lclP->V, lclP->V0);CHKERRQ(ierr);
263a7e14dcfSSatish Balay     ierr = VecCopy(lclP->GU,lclP->GU0);CHKERRQ(ierr);
264a7e14dcfSSatish Balay     ierr = VecCopy(lclP->GV,lclP->GV0);CHKERRQ(ierr);
265a7e14dcfSSatish Balay     ierr = VecCopy(lclP->GAugL_U,lclP->GAugL_U0);CHKERRQ(ierr);
266a7e14dcfSSatish Balay     ierr = VecCopy(lclP->GAugL_V,lclP->GAugL_V0);CHKERRQ(ierr);
267a7e14dcfSSatish Balay     ierr = VecCopy(lclP->GL_U,lclP->GL_U0);CHKERRQ(ierr);
268a7e14dcfSSatish Balay     ierr = VecCopy(lclP->GL_V,lclP->GL_V0);CHKERRQ(ierr);
269a7e14dcfSSatish Balay 
270a7e14dcfSSatish Balay     lclP->aug0 = lclP->aug;
271a7e14dcfSSatish Balay     lclP->lgn0 = lclP->lgn;
272a7e14dcfSSatish Balay 
273a7e14dcfSSatish Balay     /* Given the design variables, we need to project the current iterate
274a7e14dcfSSatish Balay        onto the linearized constraint.  We choose to fix the design variables
275a7e14dcfSSatish Balay        and solve the linear system for the state variables.  The resulting
276a7e14dcfSSatish Balay        point is the Newton direction */
277a7e14dcfSSatish Balay 
278a7e14dcfSSatish Balay     /* Solve r = A\con */
279a7e14dcfSSatish Balay     lclP->solve_type = LCL_FORWARD1;
280a7e14dcfSSatish Balay     ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
281a7e14dcfSSatish Balay 
282a7e14dcfSSatish Balay     if (tao->jacobian_state_inv) {
283a7e14dcfSSatish Balay       ierr = MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);CHKERRQ(ierr);
284a7e14dcfSSatish Balay     } else {
28523ee1639SBarry Smith       ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr);
286a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tao->constraints,  lclP->r);CHKERRQ(ierr);
287a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
288a7e14dcfSSatish Balay       tao->ksp_its+=its;
289ae93cb3cSJason Sarich       tao->ksp_tot_its+=tao->ksp_its;
290a7e14dcfSSatish Balay     }
291a7e14dcfSSatish Balay 
292a7e14dcfSSatish Balay     /* Set design step direction dv to zero */
293a7e14dcfSSatish Balay     ierr = VecSet(lclP->s, 0.0);CHKERRQ(ierr);
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay     /*
296a7e14dcfSSatish Balay        Check sufficient descent for constraint merit function .5*||con||^2
297a7e14dcfSSatish Balay        con' Ak r >= eps1 ||r||^(2+eps2)
298a7e14dcfSSatish Balay     */
299a7e14dcfSSatish Balay 
300a7e14dcfSSatish Balay     /* Compute WU= Ak' * con */
301a7e14dcfSSatish Balay     if (symmetric)  {
302a7e14dcfSSatish Balay       ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr);
303a7e14dcfSSatish Balay     } else {
304a7e14dcfSSatish Balay       ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr);
305a7e14dcfSSatish Balay     }
306a7e14dcfSSatish Balay     /* Compute r * Ak' * con */
307a7e14dcfSSatish Balay     ierr = VecDot(lclP->r,lclP->WU,&rWU);CHKERRQ(ierr);
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay     /* compute ||r||^(2+eps2) */
310a7e14dcfSSatish Balay     ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr);
311a7e14dcfSSatish Balay     r2 = PetscPowScalar(r2,2.0+lclP->eps2);
312a7e14dcfSSatish Balay     adec = lclP->eps1 * r2;
313a7e14dcfSSatish Balay 
314a7e14dcfSSatish Balay     if (rWU < adec) {
315955c1f14SBarry Smith       ierr = PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n");CHKERRQ(ierr);
316a7e14dcfSSatish Balay       if (lclP->verbose) {
317f06e3bfaSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);CHKERRQ(ierr);
318a7e14dcfSSatish Balay       }
319a7e14dcfSSatish Balay 
320a7e14dcfSSatish Balay       ierr = PetscInfo(tao,"Using steepest descent direction instead.\n");CHKERRQ(ierr);
321a7e14dcfSSatish Balay       ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr);
322a7e14dcfSSatish Balay       ierr = VecAXPY(lclP->r,-1.0,lclP->WU);CHKERRQ(ierr);
323a7e14dcfSSatish Balay       ierr = VecDot(lclP->r,lclP->r,&rWU);CHKERRQ(ierr);
324a7e14dcfSSatish Balay       ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr);
325a7e14dcfSSatish Balay       r2 = PetscPowScalar(r2,2.0+lclP->eps2);
326a7e14dcfSSatish Balay       ierr = VecDot(lclP->r,lclP->GAugL_U,&descent);CHKERRQ(ierr);
327a7e14dcfSSatish Balay       adec = lclP->eps1 * r2;
328a7e14dcfSSatish Balay     }
329a7e14dcfSSatish Balay 
330a7e14dcfSSatish Balay     /*
331a7e14dcfSSatish Balay        Check descent for aug. lagrangian
332a7e14dcfSSatish Balay        r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2)
333a7e14dcfSSatish Balay           GL_U = GUk - Ak'*yk
334a7e14dcfSSatish Balay           WU   = Ak'*con
335a7e14dcfSSatish Balay                                          adec=eps1||r||^(2+eps2)
336a7e14dcfSSatish Balay 
337a7e14dcfSSatish Balay        ==>
338a7e14dcfSSatish Balay        Check r'GL_U - rho*r'WU <= adec
339a7e14dcfSSatish Balay     */
340a7e14dcfSSatish Balay 
341302440fdSBarry Smith     ierr = VecDot(lclP->r,lclP->GL_U,&rGL_U);CHKERRQ(ierr);
342a7e14dcfSSatish Balay     aldescent =  rGL_U - lclP->rho*rWU;
343a7e14dcfSSatish Balay     if (aldescent > -adec) {
344a7e14dcfSSatish Balay       if (lclP->verbose) {
345f06e3bfaSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr);
346a7e14dcfSSatish Balay       }
3477d3de750SJacob Faibussowitsch       ierr = PetscInfo(tao,"Newton direction not descent for augmented Lagrangian: %g\n",(double)aldescent);CHKERRQ(ierr);
348a7e14dcfSSatish Balay       lclP->rho =  (rGL_U - adec)/rWU;
349a7e14dcfSSatish Balay       if (lclP->rho > lclP->rhomax) {
350a7e14dcfSSatish Balay         lclP->rho = lclP->rhomax;
35198921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"rho=%g > rhomax, case not implemented.  Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho);
352a7e14dcfSSatish Balay       }
353a7e14dcfSSatish Balay       if (lclP->verbose) {
354f06e3bfaSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"  Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr);
355a7e14dcfSSatish Balay       }
3567d3de750SJacob Faibussowitsch       ierr = PetscInfo(tao,"  Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr);
357a7e14dcfSSatish Balay     }
358a7e14dcfSSatish Balay 
359a7e14dcfSSatish Balay     ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
360a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
361a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
362a7e14dcfSSatish Balay 
363a7e14dcfSSatish Balay     /* We now minimize the augmented Lagrangian along the Newton direction */
364a7e14dcfSSatish Balay     ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr);
365a7e14dcfSSatish Balay     ierr = LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr);
366a7e14dcfSSatish Balay     ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr);
367a7e14dcfSSatish Balay     ierr = LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);CHKERRQ(ierr);
368a7e14dcfSSatish Balay     ierr = LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);CHKERRQ(ierr);
369a7e14dcfSSatish Balay 
370a7e14dcfSSatish Balay     lclP->recompute_jacobian_flag = PETSC_TRUE;
371a7e14dcfSSatish Balay 
372a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
3738caf6e8cSBarry Smith     ierr = TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);CHKERRQ(ierr);
374a7e14dcfSSatish Balay     ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
375a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr);
376a7e14dcfSSatish Balay     if (lclP->verbose) {
377f06e3bfaSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);CHKERRQ(ierr);
378a7e14dcfSSatish Balay     }
379a7e14dcfSSatish Balay 
380a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
381a7e14dcfSSatish Balay     ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
382a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
383a7e14dcfSSatish Balay 
384a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
385a7e14dcfSSatish Balay 
386a7e14dcfSSatish Balay     /* Check convergence */
387a7e14dcfSSatish Balay     ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
388a7e14dcfSSatish Balay     ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
389763847b4SAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
390763847b4SAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
391763847b4SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
392763847b4SAlp Dener     if (tao->reason != TAO_CONTINUE_ITERATING) {
393a7e14dcfSSatish Balay       break;
394a7e14dcfSSatish Balay     }
395a7e14dcfSSatish Balay 
396a7e14dcfSSatish Balay     /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */
397a7e14dcfSSatish Balay     for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++) {
398a7e14dcfSSatish Balay       /* We now minimize the objective function starting from the fraction of
399a7e14dcfSSatish Balay          the Newton point accepted by applying one step of a reduced-space
400a7e14dcfSSatish Balay          method.  The optimization problem is:
401a7e14dcfSSatish Balay 
402a7e14dcfSSatish Balay          minimize f(u+du, v+dv)
403a7e14dcfSSatish Balay          s. t.    A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0)
404a7e14dcfSSatish Balay 
405a7e14dcfSSatish Balay          In particular, we have that
406a7e14dcfSSatish Balay          du = -inv(A)*(Bdv + alpha g) */
407a7e14dcfSSatish Balay 
408ffad9901SBarry Smith       ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
40994ab13aaSBarry Smith       ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr);
410a7e14dcfSSatish Balay 
411a7e14dcfSSatish Balay       /* Store V and constraints */
412a7e14dcfSSatish Balay       ierr = VecCopy(lclP->V, lclP->V1);CHKERRQ(ierr);
413a7e14dcfSSatish Balay       ierr = VecCopy(tao->constraints,lclP->con1);CHKERRQ(ierr);
414a7e14dcfSSatish Balay 
415a7e14dcfSSatish Balay       /* Compute multipliers */
416a7e14dcfSSatish Balay       /* p1 */
417a7e14dcfSSatish Balay       ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
418a7e14dcfSSatish Balay       lclP->solve_type = LCL_ADJOINT1;
419a7e14dcfSSatish Balay       ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
420a7e14dcfSSatish Balay       if (tao->jacobian_state_pre) {
421a7e14dcfSSatish Balay         ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
422a7e14dcfSSatish Balay       } else {
423a7e14dcfSSatish Balay         pset = pflag = PETSC_TRUE;
424a7e14dcfSSatish Balay       }
425f06e3bfaSBarry Smith       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
426f06e3bfaSBarry Smith       else symmetric = PETSC_FALSE;
427a7e14dcfSSatish Balay 
428a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
429a7e14dcfSSatish Balay         if (symmetric) {
430a7e14dcfSSatish Balay           ierr = MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr);
431a7e14dcfSSatish Balay         } else {
432a7e14dcfSSatish Balay           ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr);
433a7e14dcfSSatish Balay         }
434a7e14dcfSSatish Balay       } else {
435a7e14dcfSSatish Balay         if (symmetric) {
436a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, lclP->GAugL_U,  lclP->lamda);CHKERRQ(ierr);
437a7e14dcfSSatish Balay         } else {
438a7e14dcfSSatish Balay           ierr = KSPSolveTranspose(tao->ksp, lclP->GAugL_U,  lclP->lamda);CHKERRQ(ierr);
439a7e14dcfSSatish Balay         }
440a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
441a7e14dcfSSatish Balay         tao->ksp_its+=its;
442ae93cb3cSJason Sarich         tao->ksp_tot_its+=its;
443a7e14dcfSSatish Balay       }
444a7e14dcfSSatish Balay       ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);CHKERRQ(ierr);
445a7e14dcfSSatish Balay       ierr = VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);CHKERRQ(ierr);
446a7e14dcfSSatish Balay 
447a7e14dcfSSatish Balay       /* Compute the limited-memory quasi-newton direction */
4488931d482SJason Sarich       if (tao->niter > 0) {
4499515a401SAlp Dener         ierr = MatSolve(lclP->R,lclP->g1,lclP->s);CHKERRQ(ierr);
450a7e14dcfSSatish Balay         ierr = VecDot(lclP->s,lclP->g1,&descent);CHKERRQ(ierr);
4510583ad50SJason Sarich         if (descent <= 0) {
452a7e14dcfSSatish Balay           if (lclP->verbose) {
453f06e3bfaSBarry Smith             ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);CHKERRQ(ierr);
454a7e14dcfSSatish Balay           }
455a7e14dcfSSatish Balay           ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr);
456a7e14dcfSSatish Balay         }
457a7e14dcfSSatish Balay       } else {
458a7e14dcfSSatish Balay         ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr);
459a7e14dcfSSatish Balay       }
460a7e14dcfSSatish Balay       ierr = VecScale(lclP->g1,-1.0);CHKERRQ(ierr);
461a7e14dcfSSatish Balay 
462a7e14dcfSSatish Balay       /* Recover the full space direction */
463a7e14dcfSSatish Balay       ierr = MatMult(tao->jacobian_design,lclP->s,lclP->WU);CHKERRQ(ierr);
464e9f9aeaeSSatish Balay       /* ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); */ /*  Initial guess in CG */
465a7e14dcfSSatish Balay       lclP->solve_type = LCL_FORWARD2;
466a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
467a7e14dcfSSatish Balay         ierr = MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);CHKERRQ(ierr);
468a7e14dcfSSatish Balay       } else {
469a7e14dcfSSatish Balay         ierr = KSPSolve(tao->ksp, lclP->WU, lclP->r);CHKERRQ(ierr);
470a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
471a7e14dcfSSatish Balay         tao->ksp_its += its;
472ae93cb3cSJason Sarich         tao->ksp_tot_its+=its;
473a7e14dcfSSatish Balay       }
474a7e14dcfSSatish Balay 
475a7e14dcfSSatish Balay       /* We now minimize the augmented Lagrangian along the direction -r,s */
476a7e14dcfSSatish Balay       ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr);
477a7e14dcfSSatish Balay       ierr = LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr);
478a7e14dcfSSatish Balay       ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr);
479a7e14dcfSSatish Balay       lclP->recompute_jacobian_flag = PETSC_TRUE;
480a7e14dcfSSatish Balay 
481a7e14dcfSSatish Balay       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
4828caf6e8cSBarry Smith       ierr = TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);CHKERRQ(ierr);
483a7e14dcfSSatish Balay       ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
484a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);CHKERRQ(ierr);
485a7e14dcfSSatish Balay       if (lclP->verbose) {
486f06e3bfaSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength =  %g\n",(double)step);CHKERRQ(ierr);
487a7e14dcfSSatish Balay       }
488a7e14dcfSSatish Balay 
489a7e14dcfSSatish Balay       ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
490a7e14dcfSSatish Balay       ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
491a7e14dcfSSatish Balay       ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
492a7e14dcfSSatish Balay       ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
493a7e14dcfSSatish Balay       ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
494a7e14dcfSSatish Balay 
495a7e14dcfSSatish Balay       /* Compute the reduced gradient at the new point */
496a7e14dcfSSatish Balay 
497ffad9901SBarry Smith       ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
49894ab13aaSBarry Smith       ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr);
499a7e14dcfSSatish Balay 
500a7e14dcfSSatish Balay       /* p2 */
501a7e14dcfSSatish Balay       /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */
502a7e14dcfSSatish Balay       if (phase2_iter==0) {
503a7e14dcfSSatish Balay         ierr = VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);CHKERRQ(ierr);
504f06e3bfaSBarry Smith       } else {
505a7e14dcfSSatish Balay         ierr = VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);CHKERRQ(ierr);
506a7e14dcfSSatish Balay       }
507a7e14dcfSSatish Balay 
508a7e14dcfSSatish Balay       ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
509a7e14dcfSSatish Balay       if (tao->jacobian_state_pre) {
510a7e14dcfSSatish Balay         ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
511a7e14dcfSSatish Balay       } else {
512a7e14dcfSSatish Balay         pset = pflag = PETSC_TRUE;
513a7e14dcfSSatish Balay       }
514f06e3bfaSBarry Smith       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
515f06e3bfaSBarry Smith       else symmetric = PETSC_FALSE;
516a7e14dcfSSatish Balay 
517a7e14dcfSSatish Balay       lclP->solve_type = LCL_ADJOINT2;
518a7e14dcfSSatish Balay       if (tao->jacobian_state_inv) {
519a7e14dcfSSatish Balay         if (symmetric) {
520a7e14dcfSSatish Balay           ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
521a7e14dcfSSatish Balay         } else {
522a7e14dcfSSatish Balay           ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
523a7e14dcfSSatish Balay         }
524a7e14dcfSSatish Balay       } else {
525a7e14dcfSSatish Balay         if (symmetric) {
526a7e14dcfSSatish Balay           ierr = KSPSolve(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
527a7e14dcfSSatish Balay         } else {
528a7e14dcfSSatish Balay           ierr = KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
529a7e14dcfSSatish Balay         }
530a7e14dcfSSatish Balay         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
531a7e14dcfSSatish Balay         tao->ksp_its += its;
5322d9aa51bSJason Sarich         tao->ksp_tot_its += its;
533a7e14dcfSSatish Balay       }
534a7e14dcfSSatish Balay 
535a7e14dcfSSatish Balay       ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);CHKERRQ(ierr);
536a7e14dcfSSatish Balay       ierr = VecAXPY(lclP->g2,-1.0,lclP->GV);CHKERRQ(ierr);
537a7e14dcfSSatish Balay 
538a7e14dcfSSatish Balay       ierr = VecScale(lclP->g2,-1.0);CHKERRQ(ierr);
539a7e14dcfSSatish Balay 
540a7e14dcfSSatish Balay       /* Update the quasi-newton approximation */
541a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(lclP->R,lclP->V,lclP->g2);CHKERRQ(ierr);
542a7e14dcfSSatish 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 */
543a7e14dcfSSatish Balay 
544a7e14dcfSSatish Balay     }
545a7e14dcfSSatish Balay 
546a7e14dcfSSatish Balay     ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr);
547a7e14dcfSSatish Balay 
548a7e14dcfSSatish Balay     /* Evaluate Function, Gradient, Constraints, and Jacobian */
549a7e14dcfSSatish Balay     ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
550a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
551a7e14dcfSSatish Balay     ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
552a7e14dcfSSatish Balay 
553ffad9901SBarry Smith     ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
55494ab13aaSBarry Smith     ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr);
555a7e14dcfSSatish Balay     ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr);
556a7e14dcfSSatish Balay 
557a7e14dcfSSatish Balay     ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
558a7e14dcfSSatish Balay 
559a7e14dcfSSatish Balay     ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
560a7e14dcfSSatish Balay 
561a7e14dcfSSatish Balay     /* Evaluate constraint norm */
562a7e14dcfSSatish Balay     ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
563a7e14dcfSSatish Balay 
564a7e14dcfSSatish Balay     /* Monitor convergence */
5658931d482SJason Sarich     tao->niter++;
566763847b4SAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
567763847b4SAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
568763847b4SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
569a7e14dcfSSatish Balay   }
570a7e14dcfSSatish Balay   PetscFunctionReturn(0);
571a7e14dcfSSatish Balay }
572a7e14dcfSSatish Balay 
5731522df2eSJason Sarich /*MC
5741522df2eSJason Sarich  TAOLCL - linearly constrained lagrangian method for pde-constrained optimization
5751522df2eSJason Sarich 
5761522df2eSJason Sarich + -tao_lcl_eps1 - epsilon 1 tolerance
57794ae4db5SBarry Smith . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr);
57894ae4db5SBarry Smith . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr);
57994ae4db5SBarry Smith . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr);
5801522df2eSJason Sarich . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm
5811522df2eSJason Sarich . -tao_lcl_verbose - Print verbose output if True
5821522df2eSJason Sarich . -tao_lcl_tola - Tolerance for first forward solve
5831522df2eSJason Sarich . -tao_lcl_tolb - Tolerance for first adjoint solve
5841522df2eSJason Sarich . -tao_lcl_tolc - Tolerance for second forward solve
5851522df2eSJason Sarich - -tao_lcl_told - Tolerance for second adjoint solve
5861522df2eSJason Sarich 
5871eb8069cSJason Sarich   Level: beginner
5881522df2eSJason Sarich M*/
589728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao)
590a7e14dcfSSatish Balay {
591a7e14dcfSSatish Balay   TAO_LCL        *lclP;
592a7e14dcfSSatish Balay   PetscErrorCode ierr;
5938caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
594a7e14dcfSSatish Balay 
595f06e3bfaSBarry Smith   PetscFunctionBegin;
596a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_LCL;
597a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_LCL;
598a7e14dcfSSatish Balay   tao->ops->view = TaoView_LCL;
599a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LCL;
600a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_LCL;
6013c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&lclP);CHKERRQ(ierr);
602a7e14dcfSSatish Balay   tao->data = (void*)lclP;
603a7e14dcfSSatish Balay 
6046552cf8aSJason Sarich   /* Override default settings (unless already changed) */
6056552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 200;
6066552cf8aSJason Sarich   if (!tao->catol_changed) tao->catol = 1.0e-4;
6076552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gttol = 1.0e-4;
6086552cf8aSJason Sarich   if (!tao->grtol_changed) tao->gttol = 1.0e-4;
6096552cf8aSJason Sarich   if (!tao->gttol_changed) tao->gttol = 1.0e-4;
610a7e14dcfSSatish Balay   lclP->rho0 = 1.0e-4;
611a7e14dcfSSatish Balay   lclP->rhomax=1e5;
612a7e14dcfSSatish Balay   lclP->eps1 = 1.0e-8;
613a7e14dcfSSatish Balay   lclP->eps2 = 0.0;
614a7e14dcfSSatish Balay   lclP->solve_type=2;
615a7e14dcfSSatish Balay   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
616a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
61763b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
618a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
6195d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
620a7e14dcfSSatish Balay 
621a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao);CHKERRQ(ierr);
622a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
62363b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
6245d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
625a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
6260c51296cSAlp Dener 
6270c51296cSAlp Dener   ierr = MatCreate(((PetscObject)tao)->comm, &lclP->R);CHKERRQ(ierr);
6280c51296cSAlp Dener   ierr = MatSetType(lclP->R, MATLMVMBFGS);CHKERRQ(ierr);
629a7e14dcfSSatish Balay   PetscFunctionReturn(0);
630a7e14dcfSSatish Balay }
631a7e14dcfSSatish Balay 
632a7e14dcfSSatish Balay static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
633a7e14dcfSSatish Balay {
634441846f8SBarry Smith   Tao            tao = (Tao)ptr;
635a7e14dcfSSatish Balay   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
636a7e14dcfSSatish Balay   PetscBool      set,pset,flag,pflag,symmetric;
637a7e14dcfSSatish Balay   PetscReal      cdotl;
638a7e14dcfSSatish Balay   PetscErrorCode ierr;
639a7e14dcfSSatish Balay 
640a7e14dcfSSatish Balay   PetscFunctionBegin;
641a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao,X,f,G);CHKERRQ(ierr);
642a7e14dcfSSatish Balay   ierr = LCLScatter(lclP,G,lclP->GU,lclP->GV);CHKERRQ(ierr);
643a7e14dcfSSatish Balay   if (lclP->recompute_jacobian_flag) {
644ffad9901SBarry Smith     ierr = TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
64594ab13aaSBarry Smith     ierr = TaoComputeJacobianDesign(tao,X,tao->jacobian_design);CHKERRQ(ierr);
646a7e14dcfSSatish Balay   }
647a7e14dcfSSatish Balay   ierr = TaoComputeConstraints(tao,X, tao->constraints);CHKERRQ(ierr);
648a7e14dcfSSatish Balay   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
649a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
650a7e14dcfSSatish Balay     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
651a7e14dcfSSatish Balay   } else {
652a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
653a7e14dcfSSatish Balay   }
654f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
655f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
656a7e14dcfSSatish Balay 
657a7e14dcfSSatish Balay   ierr = VecDot(lclP->lamda0, tao->constraints, &cdotl);CHKERRQ(ierr);
658a7e14dcfSSatish Balay   lclP->lgn = *f - cdotl;
659a7e14dcfSSatish Balay 
660a7e14dcfSSatish Balay   /* Gradient of Lagrangian GL = G - J' * lamda */
661a7e14dcfSSatish Balay   /*      WU = A' * WL
662a7e14dcfSSatish Balay           WV = B' * WL */
663a7e14dcfSSatish Balay   if (symmetric) {
664a7e14dcfSSatish Balay     ierr = MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr);
665a7e14dcfSSatish Balay   } else {
666a7e14dcfSSatish Balay     ierr = MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr);
667a7e14dcfSSatish Balay   }
668a7e14dcfSSatish Balay   ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V);CHKERRQ(ierr);
669a7e14dcfSSatish Balay   ierr = VecScale(lclP->GL_U,-1.0);CHKERRQ(ierr);
670a7e14dcfSSatish Balay   ierr = VecScale(lclP->GL_V,-1.0);CHKERRQ(ierr);
671a7e14dcfSSatish Balay   ierr = VecAXPY(lclP->GL_U,1.0,lclP->GU);CHKERRQ(ierr);
672a7e14dcfSSatish Balay   ierr = VecAXPY(lclP->GL_V,1.0,lclP->GV);CHKERRQ(ierr);
673a7e14dcfSSatish Balay   ierr = LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL);CHKERRQ(ierr);
674a7e14dcfSSatish Balay 
675a7e14dcfSSatish Balay   f[0] = lclP->lgn;
676302440fdSBarry Smith   ierr = VecCopy(lclP->GL,G);CHKERRQ(ierr);
677a7e14dcfSSatish Balay   PetscFunctionReturn(0);
678a7e14dcfSSatish Balay }
679a7e14dcfSSatish Balay 
680a7e14dcfSSatish Balay static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
681a7e14dcfSSatish Balay {
682441846f8SBarry Smith   Tao            tao = (Tao)ptr;
683a7e14dcfSSatish Balay   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
684a7e14dcfSSatish Balay   PetscReal      con2;
685a7e14dcfSSatish Balay   PetscBool      flag,pflag,set,pset,symmetric;
686a7e14dcfSSatish Balay   PetscErrorCode ierr;
687a7e14dcfSSatish Balay 
688a7e14dcfSSatish Balay   PetscFunctionBegin;
689a7e14dcfSSatish Balay   ierr = LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao);CHKERRQ(ierr);
690a7e14dcfSSatish Balay   ierr = LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
691a7e14dcfSSatish Balay   ierr = VecDot(tao->constraints,tao->constraints,&con2);CHKERRQ(ierr);
692a7e14dcfSSatish Balay   lclP->aug = lclP->lgn + 0.5*lclP->rho*con2;
693a7e14dcfSSatish Balay 
694a7e14dcfSSatish Balay   /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */
695a7e14dcfSSatish Balay   /*      WU = A' * c
696a7e14dcfSSatish Balay           WV = B' * c */
697a7e14dcfSSatish Balay   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
698a7e14dcfSSatish Balay   if (tao->jacobian_state_pre) {
699a7e14dcfSSatish Balay     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
700a7e14dcfSSatish Balay   } else {
701a7e14dcfSSatish Balay     pset = pflag = PETSC_TRUE;
702a7e14dcfSSatish Balay   }
703f06e3bfaSBarry Smith   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
704f06e3bfaSBarry Smith   else symmetric = PETSC_FALSE;
705a7e14dcfSSatish Balay 
706a7e14dcfSSatish Balay   if (symmetric) {
707a7e14dcfSSatish Balay     ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr);
708a7e14dcfSSatish Balay   } else {
709a7e14dcfSSatish Balay     ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr);
710a7e14dcfSSatish Balay   }
711a7e14dcfSSatish Balay 
712a7e14dcfSSatish Balay   ierr = MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V);CHKERRQ(ierr);
713a7e14dcfSSatish Balay   ierr = VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U);CHKERRQ(ierr);
714a7e14dcfSSatish Balay   ierr = VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V);CHKERRQ(ierr);
715a7e14dcfSSatish Balay   ierr = LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL);CHKERRQ(ierr);
716a7e14dcfSSatish Balay 
717a7e14dcfSSatish Balay   f[0] = lclP->aug;
718a7e14dcfSSatish Balay   ierr = VecCopy(lclP->GAugL,G);CHKERRQ(ierr);
719a7e14dcfSSatish Balay   PetscFunctionReturn(0);
720a7e14dcfSSatish Balay }
721a7e14dcfSSatish Balay 
722a7e14dcfSSatish Balay PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x)
723a7e14dcfSSatish Balay {
724a7e14dcfSSatish Balay   PetscErrorCode ierr;
725a7e14dcfSSatish Balay   PetscFunctionBegin;
726a7e14dcfSSatish Balay   ierr = VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
727a7e14dcfSSatish Balay   ierr = VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
728a7e14dcfSSatish Balay   ierr = VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
729a7e14dcfSSatish Balay   ierr = VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
730a7e14dcfSSatish Balay   PetscFunctionReturn(0);
731a7e14dcfSSatish Balay 
732a7e14dcfSSatish Balay }
733a7e14dcfSSatish Balay PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v)
734a7e14dcfSSatish Balay {
735a7e14dcfSSatish Balay   PetscErrorCode ierr;
736a7e14dcfSSatish Balay   PetscFunctionBegin;
737a7e14dcfSSatish Balay   ierr = VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
738a7e14dcfSSatish Balay   ierr = VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
739a7e14dcfSSatish Balay   ierr = VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
740a7e14dcfSSatish Balay   ierr = VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
741a7e14dcfSSatish Balay   PetscFunctionReturn(0);
742a7e14dcfSSatish Balay 
743a7e14dcfSSatish Balay }
744