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