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